ROMS
Loading...
Searching...
No Matches
vwalk_floats_mod Module Reference

Functions/Subroutines

subroutine, public vwalk_floats (ng, lstr, lend, predictor, my_thread, nudg)
 
subroutine vwalk_floats_tile (ng, lstr, lend, nfm3, nfm2, nfm1, nf, nfp1, predictor, my_thread, bounded, tinfo, rwalk, nudg, track)
 

Function/Subroutine Documentation

◆ vwalk_floats()

subroutine, public vwalk_floats_mod::vwalk_floats ( integer, intent(in) ng,
integer, intent(in) lstr,
integer, intent(in) lend,
logical, intent(in) predictor,
logical, dimension(lstr:), intent(in) my_thread,
real(r8), dimension(lstr:), intent(inout) nudg )

Definition at line 32 of file vwalk_floats.F.

34!***********************************************************************
35!
36 USE mod_param
37 USE mod_floats
38 USE mod_stepping
39!
40! Imported variable declarations.
41!
42 integer, intent(in) :: ng, Lstr, Lend
43
44 logical, intent(in) :: Predictor
45# ifdef ASSUMED_SHAPE
46 logical, intent(in) :: my_thread(Lstr:)
47
48 real(r8), intent(inout) :: nudg(Lstr:)
49# else
50 logical, intent(in) :: my_thread(Lstr:Lend)
51
52 real(r8), intent(inout) :: nudg(Lstr:Lend)
53# endif
54!
55! Local variable declarations.
56!
57 character (len=*), parameter :: MyFile = &
58 & __FILE__
59!
60# ifdef PROFILE
61 CALL wclock_on (ng, inlm, 10, __line__, myfile)
62# endif
63 CALL vwalk_floats_tile (ng, lstr, lend, &
64 & nfm3(ng), nfm2(ng), nfm1(ng), nf(ng), &
65 & nfp1(ng), &
66 & predictor, my_thread, &
67 & drifter(ng) % bounded, &
68 & drifter(ng) % Tinfo, &
69 & drifter(ng) % rwalk, &
70 & nudg, &
71 & drifter(ng) % track)
72# ifdef PROFILE
73 CALL wclock_off (ng, inlm, 10, __line__, myfile)
74# endif
75!
76 RETURN
type(t_drifter), dimension(:), allocatable drifter
Definition mod_floats.F:67
integer, parameter inlm
Definition mod_param.F:662
integer, dimension(:), allocatable nfm2
integer, dimension(:), allocatable nfm1
integer, dimension(:), allocatable nf
integer, dimension(:), allocatable nfm3
integer, dimension(:), allocatable nfp1
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

References mod_floats::drifter, mod_param::inlm, mod_stepping::nf, mod_stepping::nfm1, mod_stepping::nfm2, mod_stepping::nfm3, mod_stepping::nfp1, vwalk_floats_tile(), wclock_off(), and wclock_on().

Referenced by step_floats_mod::step_floats_tile().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ vwalk_floats_tile()

subroutine vwalk_floats_mod::vwalk_floats_tile ( integer, intent(in) ng,
integer, intent(in) lstr,
integer, intent(in) lend,
integer, intent(in) nfm3,
integer, intent(in) nfm2,
integer, intent(in) nfm1,
integer, intent(in) nf,
integer, intent(in) nfp1,
logical, intent(in) predictor,
logical, dimension(lstr:), intent(in) my_thread,
logical, dimension(:), intent(in) bounded,
real(r8), dimension(0:,:), intent(in) tinfo,
real(r8), dimension(:), intent(inout) rwalk,
real(r8), dimension(lstr:), intent(inout) nudg,
real(r8), dimension(:,0:,:), intent(inout) track )
private

Definition at line 81 of file vwalk_floats.F.

85!***********************************************************************
86!
87 USE mod_param
88 USE mod_parallel
89 USE mod_floats
90 USE mod_grid
91 USE mod_mixing
92 USE mod_ncparam
93 USE mod_ocean
94 USE mod_scalars
95!
97# ifdef DISTRIBUTE
98 USE distribute_mod, ONLY : mp_bcastf
99# endif
100 USE nrutil, ONLY : gasdev
101!
102! Imported variable declarations.
103!
104 integer, intent(in) :: ng, Lstr, Lend
105 integer, intent(in) :: nfm3, nfm2, nfm1, nf, nfp1
106 logical, intent(in) :: Predictor
107!
108# ifdef ASSUMED_SHAPE
109 logical, intent(in) :: bounded(:)
110 logical, intent(in) :: my_thread(Lstr:)
111
112 real(r8), intent(in) :: Tinfo(0:,:)
113
114 real(r8), intent(inout) :: rwalk(:)
115 real(r8), intent(inout) :: nudg(Lstr:)
116 real(r8), intent(inout) :: track(:,0:,:)
117# else
118 logical, intent(in) :: bounded(Nfloats(ng))
119 logical, intent(in) :: my_thread(Lstr:Lend)
120
121 real(r8), intent(in) :: Tinfo(0:izrhs,Nfloats(ng))
122
123 real(r8), intent(inout) :: rwalk(Nfloats(ng))
124 real(r8), intent(inout) :: nudg(Lstr:Lend)
125 real(r8), intent(inout) :: track(NFV(ng),0:NFT,Nfloats(ng))
126# endif
127!
128! Local variable declarations.
129!
130# ifdef MASKING
131 logical, parameter :: Lmask = .true.
132# else
133 logical, parameter :: Lmask = .false.
134# endif
135 integer :: LBi, UBi, LBj, UBj
136 integer :: i, l, nfindx
137 integer :: ierr
138
139 real(r8) :: HalfDT, akt, dakt, zrhs
140 real(r8) :: cff, cff1, cff2, cff3, cff4
141!
142! Set tile array bounds.
143!
144 lbi=lbound(grid(ng)%h,dim=1)
145 ubi=ubound(grid(ng)%h,dim=1)
146 lbj=lbound(grid(ng)%h,dim=2)
147 ubj=ubound(grid(ng)%h,dim=2)
148!
149!-----------------------------------------------------------------------
150! Compute nudging vertical velocities for random walk.
151!-----------------------------------------------------------------------
152!
153! Set float time level index to process.
154!
155 IF (predictor) THEN
156 nfindx=nf
157 ELSE
158 nfindx=nfp1
159 END IF
160!
161! If predictor step, generate random number sequence.
162!
163 IF (predictor) THEN
164# ifdef DISTRIBUTE
165 IF (master) THEN
166 CALL gasdev (rwalk)
167 END IF
168 CALL mp_bcastf (ng, inlm, rwalk)
169# else
170!$OMP MASTER
171 CALL gasdev (rwalk)
172!$OMP END MASTER
173!$OMP BARRIER
174# endif
175 END IF
176!
177! Interpolate vertical diffusion (temperature) coefficient and its
178! gradient to float locations.
179!
180 DO l=lstr,lend
181 nudg(l)=0.0_r8
182 END DO
183
184 CALL interp_floats (ng, lbi, ubi, lbj, ubj, 0, n(ng), &
185 & lstr, lend, nfindx, ifakt, isbw3d, &
186 & w3dvar, lmask, spval, nudg, &
187 & grid(ng) % pm, &
188 & grid(ng) % pn, &
189 & grid(ng) % Hz, &
190# ifdef MASKING
191 & grid(ng) % rmask, &
192# endif
193 & mixing(ng) % Akt(:,:,:,itemp), &
194 & my_thread, bounded, track)
195
196 CALL interp_floats (ng, lbi, ubi, lbj, ubj, 1, n(ng), &
197 & lstr, lend, nfindx, ifdak, isbr3d, &
198 & r3dvar, lmask, spval, nudg, &
199 & grid(ng) % pm, &
200 & grid(ng) % pn, &
201 & grid(ng) % Hz, &
202# ifdef MASKING
203 & grid(ng) % rmask, &
204# endif
205 & mixing(ng) % dAktdz, &
206 & my_thread, bounded, track)
207!
208! Compute nudging velocity coefficients. Use normally distributed
209! random numbers.
210!
211 cff=2.0_r8/dt(ng)
212 DO l=lstr,lend
213 IF (my_thread(l).and.bounded(l)) THEN
214 nudg(l)=sqrt(cff*max(0.0_r8,track(ifakt,nfindx,l)))*rwalk(l)+ &
215 & track(ifdak,nfindx,l)
216 ELSE
217 nudg(l)=0.0_r8
218 END IF
219 END DO
220!
221! Interpolate vertical slopes using nudging velocity coefficients.
222!
223 CALL interp_floats (ng, lbi, ubi, lbj, ubj, 0, n(ng), &
224 & lstr, lend, nfindx, izrhs, isbw3d, &
225 & -w3dvar, lmask, spval, nudg, &
226 & grid(ng) % pm, &
227 & grid(ng) % pn, &
228 & grid(ng) % Hz, &
229# ifdef MASKING
230 & grid(ng) % rmask, &
231# endif
232 & ocean(ng) % W, &
233 & my_thread, bounded, track)
234!
235! If newly relased float, initialize all time levels.
236!
237 halfdt=0.5_r8*dt(ng)
238
239 DO l=lstr,lend
240 IF (my_thread(l).and.bounded(l)) THEN
241 IF (time(ng)-halfdt.le.tinfo(itstr,l).and. &
242 & time(ng)+halfdt.gt.tinfo(itstr,l)) THEN
243 akt =track(ifakt,nfindx,l)
244 dakt=track(ifdak,nfindx,l)
245 zrhs=track(izrhs,nfindx,l)
246 DO i=0,nft
247 track(ifakt,i,l)=akt
248 track(ifakt,i,l)=dakt
249 track(izrhs,i,l)=zrhs
250 END DO
251 END IF
252 END IF
253 END DO
254!
255!-----------------------------------------------------------------------
256! Time step for vertical position.
257!-----------------------------------------------------------------------
258!
259! Assign predictor/corrector weights.
260!
261 IF (predictor) THEN
262 cff1=8.0_r8/3.0_r8
263 cff2=4.0_r8/3.0_r8
264 ELSE
265 cff1=9.0_r8/8.0_r8
266 cff2=1.0_r8/8.0_r8
267 cff3=3.0_r8/8.0_r8
268 cff4=6.0_r8/8.0_r8
269 END IF
270!
271! Compute new float vertical position.
272!
273# ifdef VWALK_FORWARD
274# if defined FLOAT_BIOLOGY
275 DO l=lstr,lend
276 IF (my_thread(l).and.bounded(l)) THEN
277 track(izgrd,nfp1,l)=track(izgrd,nf,l)+ &
278 & dt(ng)*(track(izrhs,nf,l)+ &
279 & track(iwbio,nf,l)* &
280 & track(i1ohz,nf,l))
281 END IF
282 END DO
283# else
284 DO l=lstr,lend
285 IF (my_thread(l).and.bounded(l)) THEN
286 track(izgrd,nfp1,l)=track(izgrd,nf,l)+ &
287 & dt(ng)*track(izrhs,nf,l)
288 END IF
289 END DO
290# endif
291# else
292# if defined FLOAT_BIOLOGY
293 IF (predictor) THEN
294 DO l=lstr,lend
295 IF (my_thread(l).and.bounded(l)) THEN
296 track(izgrd,nfp1,l)=track(izgrd,nfm3,l)+ &
297 & dt(ng)*(cff1*track(izrhs,nf ,l)- &
298 & cff2*track(izrhs,nfm1,l)+ &
299 & cff1*track(izrhs,nfm2,l)+ &
300 & cff1*track(iwbio,nf ,l)* &
301 & track(i1ohz,nf ,l)- &
302 & cff2*track(iwbio,nfm1,l)* &
303 & track(i1ohz,nfm1,l)+ &
304 & cff1*track(iwbio,nfm2,l)* &
305 & track(i1ohz,nfm2,l))
306 END IF
307 END DO
308 ELSE
309 DO l=lstr,lend
310 IF (my_thread(l).and.bounded(l)) THEN
311 track(izgrd,nfp1,l)=cff1*track(izgrd,nf ,l)- &
312 & cff2*track(izgrd,nfm2,l)+ &
313 & dt(ng)*(cff3*track(izrhs,nfp1,l)+ &
314 & cff4*track(izrhs,nf ,l)- &
315 & cff3*track(izrhs,nfm1,l)+ &
316 & cff3*track(iwbio,nfp1,l)* &
317 & track(i1ohz,nfp1,l)+ &
318 & cff4*track(iwbio,nf ,l)* &
319 & track(i1ohz,nf ,l)- &
320 & cff3*track(iwbio,nfm1,l)* &
321 & track(i1ohz,nf ,l))
322 END IF
323 END DO
324 END IF
325# else
326 IF (predictor) THEN
327 DO l=lstr,lend
328 IF (my_thread(l).and.bounded(l)) THEN
329 track(izgrd,nfp1,l)=track(izgrd,nfm3,l)+ &
330 & dt(ng)*(cff1*track(izrhs,nf ,l)- &
331 & cff2*track(izrhs,nfm1,l)+ &
332 & cff1*track(izrhs,nfm2,l))
333 END IF
334 END DO
335 ELSE
336 DO l=lstr,lend
337 IF (my_thread(l).and.bounded(l)) THEN
338 track(izgrd,nfp1,l)=cff1*track(izgrd,nf ,l)- &
339 & cff2*track(izgrd,nfm2,l)+ &
340 & dt(ng)*(cff3*track(izrhs,nfp1,l)+ &
341 & cff4*track(izrhs,nf ,l)- &
342 & cff3*track(izrhs,nfm1,l))
343 END IF
344 END DO
345 END IF
346# endif
347# endif
348!
349! Zeroth-out nudging velocities coefficients.
350!
351 DO l=lstr,lend
352 nudg(l)=0.0_r8
353 END DO
354!
355 RETURN
subroutine, public interp_floats(ng, lbi, ubi, lbj, ubj, lbk, ubk, lstr, lend, itime, ifield, isbval, gtype, maskit, fspval, nudg, pm, pn, hz, amask, a, my_thread, bounded, track)
integer, parameter i1ohz
Definition mod_floats.F:97
integer, parameter iwbio
Definition mod_floats.F:101
integer, parameter itstr
Definition mod_floats.F:80
integer, parameter izrhs
Definition mod_floats.F:89
integer, parameter izgrd
Definition mod_floats.F:83
integer, parameter ifakt
Definition mod_floats.F:92
integer, parameter ifdak
Definition mod_floats.F:93
type(t_grid), dimension(:), allocatable grid
Definition mod_grid.F:365
type(t_mixing), dimension(:), allocatable mixing
Definition mod_mixing.F:399
integer isbw3d
integer isbr3d
type(t_ocean), dimension(:), allocatable ocean
Definition mod_ocean.F:351
logical master
integer, dimension(:), allocatable n
Definition mod_param.F:479
integer, parameter r3dvar
Definition mod_param.F:721
integer, parameter nft
Definition mod_param.F:539
integer, parameter w3dvar
Definition mod_param.F:724
real(dp), parameter spval
real(dp), dimension(:), allocatable dt
integer itemp
real(dp), dimension(:), allocatable time
Definition nrutil.F:1

References mod_scalars::dt, mod_grid::grid, mod_floats::i1ohz, mod_floats::ifakt, mod_floats::ifdak, mod_param::inlm, interp_floats_mod::interp_floats(), mod_ncparam::isbr3d, mod_ncparam::isbw3d, mod_scalars::itemp, mod_floats::itstr, mod_floats::iwbio, mod_floats::izgrd, mod_parallel::master, mod_mixing::mixing, mod_param::n, mod_ocean::ocean, mod_param::r3dvar, mod_scalars::spval, mod_scalars::time, and mod_param::w3dvar.

Referenced by vwalk_floats().

Here is the call graph for this function:
Here is the caller graph for this function: