ROMS
Loading...
Searching...
No Matches
tl_get_idata.F
Go to the documentation of this file.
1#include "cppdefs.h"
2#ifdef TANGENT
3 SUBROUTINE tl_get_idata (ng)
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! This routine reads input data that needs to be obtained only once. !
13! !
14! Currently, this routine is only executed in serial mode by the !
15! main thread. !
16! !
17!=======================================================================
18!
19 USE mod_param
20 USE mod_grid
21 USE mod_iounits
22 USE mod_ncparam
23 USE mod_parallel
24 USE mod_scalars
25 USE mod_sources
26 USE mod_stepping
27# if defined SSH_TIDES_NOT_YET || defined UV_TIDES_NOT_YET
28 USE mod_tides
29# endif
30!
31 USE strings_mod, ONLY : founderror
32!
33 implicit none
34!
35! Imported variable declarations.
36!
37 integer, intent(in) :: ng
38!
39! Local variable declarations.
40!
41 logical, save :: recordless = .true.
42
43 logical, dimension(3) :: update = &
44 & (/ .FALSE., .FALSE., .FALSE. /)
45!
46 integer :: LBi, UBi, LBj, UBj
47 integer :: itrc, is
48!
49 real(r8) :: time_save = 0.0_r8
50!
51 character (len=*), parameter :: MyFile = &
52 & __FILE__
53!
54 sourcefile=myfile
55!
56! Lower and upper bounds for tiled arrays.
57!
58 lbi=lbound(grid(ng)%h,dim=1)
59 ubi=ubound(grid(ng)%h,dim=1)
60 lbj=lbound(grid(ng)%h,dim=2)
61 ubj=ubound(grid(ng)%h,dim=2)
62
63# ifdef PROFILE
64!
65!-----------------------------------------------------------------------
66! Turn on input data time wall clock.
67!-----------------------------------------------------------------------
68!
69 CALL wclock_on (ng, itlm, 3, __line__, myfile)
70# endif
71
72# if defined SSH_TIDES_NOT_YET || defined UV_TIDES_NOT_YET
73!
74!-----------------------------------------------------------------------
75! Tide period, amplitude, phase, and currents.
76!-----------------------------------------------------------------------
77!
78! Tidal Period.
79!
80 IF (lprocesstides(ng)) THEN
81 IF (iic(ng).eq.0) THEN
82 CALL get_ngfld (ng, itlm, idtper, tide(ng)%ncid, &
83# if defined PIO_LIB && defined DISTRIBUTE
84 & tide(ng)%pioFile, &
85# endif
86 & 1, tide(ng), recordless, update(1), &
87 & 1, mtc, 1, 1, 1, ntc(ng), 1, &
88 & tides(ng) % Tperiod)
89 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
90 END IF
91 END IF
92# endif
93
94# ifdef SSH_TIDES_NOT_YET
95!
96! Tidal elevation amplitude and phase. In order to read data as a
97! function of tidal period, we need to reset the model time variables
98! temporarily.
99!
100 IF (lprocesstides(ng)) THEN
101 IF (iic(ng).eq.0) THEN
102 time_save=time(ng)
103 time(ng)=8640000.0_r8
104 tdays(ng)=time(ng)*sec2day
105!
106 CALL get_2dfld (ng, itlm, idtzam, tide(ng)%ncid, &
107# if defined PIO_LIB && defined DISTRIBUTE
108 & tide(ng)%pioFile, &
109# endif
110 & 1, tide(ng), update(1), &
111 & lbi, ubi, lbj, ubj, mtc, ntc(ng), &
112# ifdef MASKING
113 & grid(ng) % rmask, &
114# endif
115 & tides(ng) % SSH_Tamp)
116 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
117!
118 CALL get_2dfld (ng, itlm, idtzph, tide(ng)%ncid, &
119# if defined PIO_LIB && defined DISTRIBUTE
120 & tide(ng)%pioFile, &
121# endif
122 & 1, tide(ng), update(1), &
123 & lbi, ubi, lbj, ubj, mtc, ntc(ng), &
124# ifdef MASKING
125 & grid(ng) % rmask, &
126# endif
127 & tides(ng) % SSH_Tphase)
128 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
129!
130 time(ng)=time_save
131 tdays(ng)=time(ng)*sec2day
132 END IF
133 END IF
134# endif
135
136# ifdef UV_TIDES_NOT_YET
137!
138! Tidal currents angle, phase, major and minor ellipse axis.
139!
140 IF (lprocesstides(ng)) THEN
141 IF (iic(ng).eq.0) THEN
142 time_save=time(ng)
143 time(ng)=8640000.0_r8
144 tdays(ng)=time(ng)*sec2day
145!
146 CALL get_2dfld (ng, itlm, idtvan, tide(ng)%ncid, &
147# if defined PIO_LIB && defined DISTRIBUTE
148 & tide(ng)%pioFile, &
149# endif
150 & 1, tide(ng), update(1), &
151 & lbi, ubi, lbj, ubj, mtc, ntc(ng), &
152# ifdef MASKING
153 & grid(ng) % rmask, &
154# endif
155 & tides(ng) % UV_Tangle)
156 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
157!
158 CALL get_2dfld (ng, itlm, idtvph, tide(ng)%ncid, &
159# if defined PIO_LIB && defined DISTRIBUTE
160 & tide(ng)%pioFile, &
161# endif
162 & 1, tide(ng), update(1), &
163 & lbi, ubi, lbj, ubj, mtc, ntc(ng), &
164# ifdef MASKING
165 & grid(ng) % rmask, &
166# endif
167 & tides(ng) % UV_Tphase)
168 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
169!
170 CALL get_2dfld (ng, itlm, idtvma, tide(ng)%ncid, &
171# if defined PIO_LIB && defined DISTRIBUTE
172 & tide(ng)%pioFile, &
173# endif
174 & 1, tide(ng), update(1), &
175 & lbi, ubi, lbj, ubj, mtc, ntc(ng), &
176# ifdef MASKING
177 & grid(ng) % rmask, &
178# endif
179 & tides(ng) % UV_Tmajor)
180 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
181!
182 CALL get_2dfld (ng, itlm, idtvmi, tide(ng)%ncid, &
183# if defined PIO_LIB && defined DISTRIBUTE
184 & tide(ng)%pioFile, &
185# endif
186 & 1, tide(ng), update(1), &
187 & lbi, ubi, lbj, ubj, mtc, ntc(ng), &
188# ifdef MASKING
189 & grid(ng) % rmask, &
190# endif
191 & tides(ng) % UV_Tminor)
192 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
193!
194 time(ng)=time_save
195 tdays(ng)=time(ng)*sec2day
196 END IF
197 END IF
198# endif
199
200# ifndef ANA_PSOURCE
201!
202!-----------------------------------------------------------------------
203! Read in point Sources/Sinks position, direction, special flag, and
204! mass transport nondimensional shape profile. Point sources are at
205! U- and V-points.
206!-----------------------------------------------------------------------
207!
208 IF ((iic(ng).eq.0).and. &
209 & (luvsrc(ng).or.lwsrc(ng).or.any(ltracersrc(:,ng)))) THEN
210 CALL get_ngfld (ng, itlm, idrxpo, ssf(ng)%ncid, &
211# if defined PIO_LIB && defined DISTRIBUTE
212 & ssf(ng)%pioFile, &
213# endif
214 & 1, ssf(ng), recordless, update(1), &
215 & 1, nsrc(ng), 1, 1, 1, nsrc(ng), 1, &
216 & sources(ng) % Xsrc)
217 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
218!
219 CALL get_ngfld (ng, itlm, idrepo, ssf(ng)%ncid, &
220# if defined PIO_LIB && defined DISTRIBUTE
221 & ssf(ng)%pioFile, &
222# endif
223 & 1, ssf(ng), recordless, update(1), &
224 & 1, nsrc(ng), 1, 1, 1, nsrc(ng), 1, &
225 & sources(ng) % Ysrc)
226 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
227!
228 CALL get_ngfld (ng, itlm, idrdir, ssf(ng)%ncid, &
229# if defined PIO_LIB && defined DISTRIBUTE
230 & ssf(ng)%pioFile, &
231# endif
232 & 1, ssf(ng), recordless, update(1), &
233 & 1, nsrc(ng), 1, 1, 1, nsrc(ng), 1, &
234 & sources(ng) % Dsrc)
235 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
236
237# ifdef SOLVE3D
238!
239 CALL get_ngfld (ng, itlm, idrvsh, ssf(ng)%ncid, &
240# if defined PIO_LIB && defined DISTRIBUTE
241 & ssf(ng)%pioFile, &
242# endif
243 & 1, ssf(ng), recordless, update(1), &
244 & 1, nsrc(ng), n(ng), 1, 1, nsrc(ng), n(ng), &
245 & sources(ng) % Qshape)
246 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
247# endif
248!
249 DO is=1,nsrc(ng)
250 sources(ng)%Isrc(is)= &
251 & max(1,min(nint(sources(ng)%Xsrc(is)),lm(ng)+1))
252 sources(ng)%Jsrc(is)= &
253 & max(1,min(nint(sources(ng)%Ysrc(is)),mm(ng)+1))
254 END DO
255 END IF
256# endif
257
258# ifdef PROFILE
259!
260!-----------------------------------------------------------------------
261! Turn off input data time wall clock.
262!-----------------------------------------------------------------------
263!
264 CALL wclock_off (ng, itlm, 3, __line__, myfile)
265# endif
266!
267 RETURN
268 END SUBROUTINE tl_get_idata
269#else
270 SUBROUTINE tl_get_idata
271 RETURN
272 END SUBROUTINE tl_get_idata
273#endif
subroutine get_2dfld(ng, model, ifield, ncid, piofile, nfiles, s, update, lbi, ubi, lbj, ubj, iout, irec, fmask, fout)
Definition get_2dfld.F:12
subroutine get_ngfld(ng, model, ifield, ncid, piofile, nfiles, s, recordless, update, lbi, ubi, ubj, ubk, istr, iend, jrec, fout)
Definition get_ngfld.F:9
type(t_grid), dimension(:), allocatable grid
Definition mod_grid.F:365
type(t_io), dimension(:), allocatable ssf
type(t_io), dimension(:), allocatable tide
character(len=256) sourcefile
integer idtzph
integer idtper
integer idtvmi
integer idrepo
integer idtzam
integer idtvph
integer idtvan
integer idrvsh
integer idrxpo
integer idtvma
integer idrdir
integer, dimension(:), allocatable n
Definition mod_param.F:479
integer mtc
Definition mod_param.F:564
integer, dimension(:), allocatable lm
Definition mod_param.F:455
integer, parameter itlm
Definition mod_param.F:663
integer, dimension(:), allocatable mm
Definition mod_param.F:456
logical, dimension(:), allocatable luvsrc
logical, dimension(:,:), allocatable ltracersrc
integer, dimension(:), allocatable iic
logical, dimension(:), allocatable lprocesstides
real(dp), dimension(:), allocatable tdays
logical, dimension(:), allocatable lwsrc
real(dp), parameter sec2day
integer exit_flag
real(dp), dimension(:), allocatable time
integer noerror
type(t_sources), dimension(:), allocatable sources
Definition mod_sources.F:90
integer, dimension(:), allocatable nsrc
Definition mod_sources.F:97
integer, dimension(:), allocatable ntc
type(t_tides), dimension(:), allocatable tides
Definition mod_tides.F:133
logical function, public founderror(flag, noerr, line, routine)
Definition strings.F:52
recursive subroutine wclock_off(ng, model, region, line, routine)
Definition timers.F:148
recursive subroutine wclock_on(ng, model, region, line, routine)
Definition timers.F:3
subroutine tl_get_idata(ng)
Definition tl_get_idata.F:4