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

Functions/Subroutines

subroutine, public step_floats (ng, lstr, lend)
 
subroutine step_floats_tile (ng, lstr, lend, knew, nnew, nfm3, nfm2, nfm1, nf, nfp1, bounded, ftype, tinfo, fz0, stuck, track)
 

Function/Subroutine Documentation

◆ step_floats()

subroutine, public step_floats_mod::step_floats ( integer, intent(in) ng,
integer, intent(in) lstr,
integer, intent(in) lend )

Definition at line 42 of file step_floats.F.

43!***********************************************************************
44!
45 USE mod_param
46 USE mod_floats
47 USE mod_stepping
48!
49! Imported variable declarations.
50!
51 integer, intent(in) :: ng, Lstr, Lend
52!
53! Local variable declarations.
54!
55 character (len=*), parameter :: MyFile = &
56 & __FILE__
57!
58# ifdef PROFILE
59 CALL wclock_on (ng, inlm, 10, __line__, myfile)
60# endif
61 CALL step_floats_tile (ng, lstr, lend, &
62 & knew(ng), nnew(ng), nfm3(ng), nfm2(ng), &
63 & nfm1(ng), nf(ng), nfp1(ng), &
64 & drifter(ng) % bounded, &
65 & drifter(ng) % Ftype, &
66 & drifter(ng) % Tinfo, &
67 & drifter(ng) % Fz0, &
68# if defined SOLVE3D && defined FLOAT_STICKY
69 & drifter(ng) % stuck, &
70# endif
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 knew
integer, dimension(:), allocatable nfm1
integer, dimension(:), allocatable nf
integer, dimension(:), allocatable nfm3
integer, dimension(:), allocatable nnew
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::knew, mod_stepping::nf, mod_stepping::nfm1, mod_stepping::nfm2, mod_stepping::nfm3, mod_stepping::nfp1, mod_stepping::nnew, step_floats_tile(), wclock_off(), and wclock_on().

Referenced by main3d().

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

◆ step_floats_tile()

subroutine step_floats_mod::step_floats_tile ( integer, intent(in) ng,
integer, intent(in) lstr,
integer, intent(in) lend,
integer, intent(in) knew,
integer, intent(in) nnew,
integer, intent(in) nfm3,
integer, intent(in) nfm2,
integer, intent(in) nfm1,
integer, intent(in) nf,
integer, intent(in) nfp1,
logical, dimension(:), intent(inout) bounded,
integer, dimension(:), intent(in) ftype,
real(r8), dimension(0:,:), intent(in) tinfo,
real(r8), dimension(:), intent(in) fz0,
logical, dimension(:), intent(inout) stuck,
real(r8), dimension(:,0:,:), intent(inout) track )
private

Definition at line 80 of file step_floats.F.

88!***********************************************************************
89!
90 USE mod_param
91 USE mod_parallel
92 USE mod_floats
93 USE mod_grid
94 USE mod_iounits
95 USE mod_ncparam
96 USE mod_ocean
97 USE mod_scalars
98!
99# ifdef FLOAT_BIOLOGY
101# endif
102# ifdef DISTRIBUTE
103 USE distribute_mod, ONLY : mp_collect
104# endif
106# if defined SOLVE3D && defined FLOAT_VWALK
108# endif
109!
110! Imported variable declarations.
111!
112 integer, intent(in) :: ng, Lstr, Lend
113 integer, intent(in) :: knew, nnew, nfm3, nfm2, nfm1, nf, nfp1
114!
115# ifdef ASSUMED_SHAPE
116 integer, intent(in) :: Ftype(:)
117 real(r8), intent(in) :: Tinfo(0:,:)
118 real(r8), intent(in) :: Fz0(:)
119
120 logical, intent(inout) :: bounded(:)
121# if defined SOLVE3D && defined FLOAT_STICKY
122 logical, intent(inout) :: stuck(:)
123# endif
124 real(r8), intent(inout) :: track(:,0:,:)
125# else
126 integer, intent(in) :: Ftype(Nfloats(ng))
127 real(r8), intent(in) :: Tinfo(0:izrhs,Nfloats(ng))
128 real(r8), intent(in) :: Fz0(Nfloats(ng))
129
130 logical, intent(inout) :: bounded(Nfloats(ng))
131# if defined SOLVE3D && defined FLOAT_STICKY
132 logical, intent(inout) :: stuck(Nfloats(ng))
133# endif
134 real(r8), intent(inout) :: track(NFV(ng),0:NFT,Nfloats(ng))
135# endif
136!
137! Local variable declarations.
138!
139 logical, parameter :: Gmask = .false.
140# ifdef MASKING
141 logical, parameter :: Lmask = .true.
142# else
143 logical, parameter :: Lmask = .false.
144# endif
145 logical, dimension(Lstr:Lend) :: my_thread
146
147 integer :: LBi, UBi, LBj, UBj
148 integer :: Ir, Jr, Npts, i, i1, i2, j, j1, j2, itrc, l, k
149
150 real(r8), parameter :: Fspv = 0.0_r8
151
152 real(r8) :: cff1, cff2, cff3, cff4, cff5, cff6, cff7, cff8, cff9
153 real(r8) :: oHz, p1, p2, q1, q2, xrhs, yrhs, zrhs, zfloat
154 real(r8) :: HalfDT
155
156 real(r8), dimension(Lstr:Lend) :: nudg
157
158# ifdef DISTRIBUTE
159 real(r8) :: Xstr, Xend, Ystr, Yend
160 real(r8), dimension(Nfloats(ng)*NFV(ng)*(NFT+1)) :: Fwrk
161# endif
162!
163! Set tile array bounds.
164!
165 lbi=lbound(grid(ng)%h,dim=1)
166 ubi=ubound(grid(ng)%h,dim=1)
167 lbj=lbound(grid(ng)%h,dim=2)
168 ubj=ubound(grid(ng)%h,dim=2)
169
170# ifdef DISTRIBUTE
171!
172!-----------------------------------------------------------------------
173! In distributed-memory configuration, determine which node bounds the
174! current location of the floats. Assign unbounded floats to the master
175! node.
176!-----------------------------------------------------------------------
177!
178! The strategy here is to build a switch that processes only the floats
179! contained within the tile node. The trajectory data for unbounded
180! floats is initialized to Fspv. These values are used during the
181! collection step at the end of the routine. Since a SUM reduction is
182! carried-out, setting Fspv to zero means the floats only contribute in
183! their own tile.
184!
185 npts=nfv(ng)*(nft+1)*nfloats(ng)
186
187 xstr=real(bounds(ng)%Istr(myrank),r8)-0.5_r8
188 xend=real(bounds(ng)%Iend(myrank),r8)+0.5_r8
189 ystr=real(bounds(ng)%Jstr(myrank),r8)-0.5_r8
190 yend=real(bounds(ng)%Jend(myrank),r8)+0.5_r8
191 DO l=lstr,lend
192 my_thread(l)=.false.
193 IF ((xstr.le.track(ixgrd,nf,l)).and. &
194 & (track(ixgrd,nf,l).lt.xend).and. &
195 & (ystr.le.track(iygrd,nf,l)).and. &
196 & (track(iygrd,nf,l).lt.yend)) THEN
197 my_thread(l)=.true.
198 ELSE IF (master.and.(.not.bounded(l))) THEN
199 my_thread(l)=.true.
200 ELSE
201 DO j=0,nft
202 DO i=1,nfv(ng)
203 track(i,j,l)=fspv
204 END DO
205 END DO
206 END IF
207 END DO
208# else
209!
210!-----------------------------------------------------------------------
211! Initialize.
212!-----------------------------------------------------------------------
213!
214 DO l=lstr,lend
215 my_thread(l)=.true.
216 END DO
217# endif
218# if !(defined SOLVE3D && defined FLOAT_VWALK)
219 DO l=lstr,lend
220 nudg(l)=0.0_r8
221 END DO
222# endif
223# if defined SOLVE3D && defined FLOAT_VWALK
224!
225!-----------------------------------------------------------------------
226! Compute vertical positions due to vertical random walk, predictor
227! step.
228!-----------------------------------------------------------------------
229!
230 CALL vwalk_floats (ng, lstr, lend, .true., my_thread, nudg)
231# endif
232!
233!-----------------------------------------------------------------------
234! Predictor step: compute first guess floats locations using a
235! 4th-order Milne time-stepping scheme.
236!-----------------------------------------------------------------------
237!
238 cff1=8.0_r8/3.0_r8
239 cff2=4.0_r8/3.0_r8
240 DO l=lstr,lend
241 IF (my_thread(l).and.bounded(l)) THEN
242 track(ixgrd,nfp1,l)=track(ixgrd,nfm3,l)+ &
243 & dt(ng)*(cff1*track(ixrhs,nf ,l)- &
244 & cff2*track(ixrhs,nfm1,l)+ &
245 & cff1*track(ixrhs,nfm2,l))
246 track(iygrd,nfp1,l)=track(iygrd,nfm3,l)+ &
247 & dt(ng)*(cff1*track(iyrhs,nf ,l)- &
248 & cff2*track(iyrhs,nfm1,l)+ &
249 & cff1*track(iyrhs,nfm2,l))
250
251# if defined SOLVE3D && !defined FLOAT_VWALK
252!
253! Compute vertical position (grid units) 3D Lagrangian floats.
254!
255 IF (ftype(l).eq.flt_lagran) THEN
256# ifdef FLOAT_BIOLOGY
257 track(izgrd,nfp1,l)=track(izgrd,nfm3,l)+ &
258 & dt(ng)*(cff1*track(izrhs,nf ,l)- &
259 & cff2*track(izrhs,nfm1,l)+ &
260 & cff1*track(izrhs,nfm2,l)+ &
261 & cff1*track(iwbio,nf ,l)* &
262 & track(i1ohz,nf ,l)- &
263 & cff2*track(iwbio,nfm1,l)* &
264 & track(i1ohz,nfm1,l)+ &
265 & cff1*track(iwbio,nfm2,l)* &
266 & track(i1ohz,nfm2,l))
267# else
268 track(izgrd,nfp1,l)=track(izgrd,nfm3,l)+ &
269 & dt(ng)*(cff1*track(izrhs,nf ,l)- &
270 & cff2*track(izrhs,nfm1,l)+ &
271 & cff1*track(izrhs,nfm2,l))
272# endif
273!
274! Compute vertical position (grid units) for isobaric floats
275! (p=g*(z+zeta)=constant) or geopotential floats (constant depth).
276! Use bilinear interpolation to determine vertical position.
277!
278 ELSE IF ((ftype(l).eq.flt_isobar).or. &
279 & (ftype(l).eq.flt_geopot)) THEN
280 ir=int(track(ixgrd,nfp1,l))
281 jr=int(track(iygrd,nfp1,l))
282!
283 i1=min(max(ir ,0),lm(ng)+1)
284 i2=min(max(ir+1,1),lm(ng)+1)
285 j1=min(max(jr ,0),mm(ng)+1)
286 j2=min(max(jr+1,0),mm(ng)+1)
287!
288 p2=real(i2-i1,r8)*(track(ixgrd,nfp1,l)-real(i1,r8))
289 q2=real(j2-j1,r8)*(track(iygrd,nfp1,l)-real(j1,r8))
290 p1=1.0_r8-p2
291 q1=1.0_r8-q2
292# ifdef MASKING
293 cff7=p1*q1*grid(ng)%z_w(i1,j1,n(ng))*grid(ng)%rmask(i1,j1)+ &
294 & p2*q1*grid(ng)%z_w(i2,j1,n(ng))*grid(ng)%rmask(i2,j1)+ &
295 & p1*q2*grid(ng)%z_w(i1,j2,n(ng))*grid(ng)%rmask(i1,j2)+ &
296 & p2*q2*grid(ng)%z_w(i2,j2,n(ng))*grid(ng)%rmask(i2,j2)
297 cff8=p1*q1*grid(ng)%rmask(i1,j1)+ &
298 & p2*q1*grid(ng)%rmask(i2,j1)+ &
299 & p1*q2*grid(ng)%rmask(i1,j2)+ &
300 & p2*q2*grid(ng)%rmask(i2,j2)
301 cff9=0.0_r8
302 IF (cff8.gt.0.0_r8) cff9=cff7/cff8
303# else
304 cff9=p1*q1*grid(ng)%z_w(i1,j1,n(ng))+ &
305 & p2*q1*grid(ng)%z_w(i2,j1,n(ng))+ &
306 & p1*q2*grid(ng)%z_w(i1,j2,n(ng))+ &
307 & p2*q2*grid(ng)%z_w(i2,j2,n(ng))
308# endif
309 cff6=cff9
310!
311 IF (ftype(l).eq.flt_geopot) THEN
312 zfloat=fz0(l)
313 ELSE IF (ftype(l).eq.flt_isobar) THEN
314 zfloat=fz0(l)+cff9
315 END IF
316!
317 DO k=n(ng)-1,0,-1
318# ifdef MASKING
319 cff7=p1*q1*grid(ng)%z_w(i1,j1,k)*grid(ng)%rmask(i1,j1)+ &
320 & p2*q1*grid(ng)%z_w(i2,j1,k)*grid(ng)%rmask(i2,j1)+ &
321 & p1*q2*grid(ng)%z_w(i1,j2,k)*grid(ng)%rmask(i1,j2)+ &
322 & p2*q2*grid(ng)%z_w(i2,j2,k)*grid(ng)%rmask(i2,j2)
323 cff8=p1*q1*grid(ng)%rmask(i1,j1)+ &
324 & p2*q1*grid(ng)%rmask(i2,j1)+ &
325 & p1*q2*grid(ng)%rmask(i1,j2)+ &
326 & p2*q2*grid(ng)%rmask(i2,j2)
327 IF (cff8.gt.0.0_r8) THEN
328 cff5=cff7/cff8
329 ELSE
330 cff5=0.0_r8
331 END IF
332# else
333 cff5=p1*q1*grid(ng)%z_w(i1,j1,k)+ &
334 & p2*q1*grid(ng)%z_w(i2,j1,k)+ &
335 & p1*q2*grid(ng)%z_w(i1,j2,k)+ &
336 & p2*q2*grid(ng)%z_w(i2,j2,k)
337# endif
338 IF ((zfloat-cff5)*(cff6-zfloat).ge.0.0_r8) THEN
339 track(izgrd,nfp1,l)=real(k,r8)+(zfloat-cff5)/(cff6-cff5)
340 END IF
341 cff6=cff5
342 END DO
343 END IF
344# endif
345 END IF
346 END DO
347!
348!-----------------------------------------------------------------------
349! Calculate slopes at new time-step.
350!-----------------------------------------------------------------------
351!
352# ifdef SOLVE3D
353 CALL interp_floats (ng, lbi, ubi, lbj, ubj, 1, n(ng), &
354 & lstr, lend, nfp1, ixrhs, isuvel, &
355 & -u3dvar, lmask, spval, nudg, &
356 & grid(ng) % pm, &
357 & grid(ng) % pn, &
358 & grid(ng) % Hz, &
359# ifdef MASKING
360 & grid(ng) % rmask, &
361# endif
362 & ocean(ng) % u(:,:,:,nnew), &
363 & my_thread, bounded, track)
364
365 CALL interp_floats (ng, lbi, ubi, lbj, ubj, 1, n(ng), &
366 & lstr, lend, nfp1, iyrhs, isvvel, &
367 & -v3dvar, lmask, spval, nudg, &
368 & grid(ng) % pm, &
369 & grid(ng) % pn, &
370 & grid(ng) % Hz, &
371# ifdef MASKING
372 & grid(ng) % rmask, &
373# endif
374 & ocean(ng) % v(:,:,:,nnew), &
375 & my_thread, bounded, track)
376
377# if !defined FLOAT_VWALK
378 CALL interp_floats (ng, lbi, ubi, lbj, ubj, 0, n(ng), &
379 & lstr, lend, nfp1, izrhs, isbw3d, &
380 & -w3dvar, lmask, spval, nudg, &
381 & grid(ng) % pm, &
382 & grid(ng) % pn, &
383 & grid(ng) % Hz, &
384# ifdef MASKING
385 & grid(ng) % rmask, &
386# endif
387 & ocean(ng) % W, &
388 & my_thread, bounded, track)
389# endif
390# if defined FLOAT_BIOLOGY
391 CALL interp_floats (ng, lbi, ubi, lbj, ubj, 1, n(ng), &
392 & lstr, lend, nfp1, i1ohz, isbr3d, &
393 & r3dvar, gmask, fspv, nudg, &
394 & grid(ng) % pm, &
395 & grid(ng) % pn, &
396 & grid(ng) % Hz, &
397# ifdef MASKING
398 & grid(ng) % rmask, &
399# endif
400 & grid(ng) % Hz, &
401 & my_thread, bounded, track)
402 DO l=lstr,lend
403 IF (my_thread(l).and.bounded(l)) THEN
404 track(i1ohz,nfp1,l)=1.0_r8/track(i1ohz,nfp1,l)
405 END IF
406 END DO
407# endif
408# else
409 CALL interp_floats (ng, lbi, ubi, lbj, ubj, 1, 1, &
410 & lstr, lend, nfp1, ixrhs, isubar, &
411 & -u2dvar, lmask, spval, nudg, &
412 & grid(ng) % pm, &
413 & grid(ng) % pn, &
414# ifdef MASKING
415 & grid(ng) % rmask, &
416# endif
417 & ocean(ng) % ubar(:,:,knew), &
418 & my_thread, bounded, track)
419
420 CALL interp_floats (ng, lbi, ubi, lbj, ubj, 1, 1, &
421 & lstr, lend, nfp1, iyrhs, isvbar, &
422 & -v2dvar, lmask, spval, nudg, &
423 & grid(ng) % pm, &
424 & grid(ng) % pn, &
425# ifdef MASKING
426 & grid(ng) % rmask, &
427# endif
428 & ocean(ng) % vbar(:,:,knew), &
429 & my_thread, bounded, track)
430# endif
431
432# ifdef FLOAT_BIOLOGY
433!
434!-----------------------------------------------------------------------
435! Compute biological float behavior, predictor step.
436!-----------------------------------------------------------------------
437!
438! Interpolate tracer to the predictor step locations. These values
439! are used in the "biology_floats" routine.
440!
441 DO itrc=1,nt(ng)
442 CALL interp_floats (ng, lbi, ubi, lbj, ubj, 1, n(ng), &
443 & lstr, lend, nfp1, iftvar(itrc), &
444 & istvar(itrc), r3dvar, lmask, spval, nudg, &
445 & grid(ng) % pm, &
446 & grid(ng) % pn, &
447 & grid(ng) % Hz, &
448# ifdef MASKING
449 & grid(ng) % rmask, &
450# endif
451 & ocean(ng) % t(:,:,:,nnew,itrc), &
452 & my_thread, bounded, track)
453 END DO
454!
455! Biological behavior predictor step.
456!
457 CALL biology_floats (ng, lstr, lend, .true., my_thread)
458# endif
459!
460!-----------------------------------------------------------------------
461! Corrector step: correct floats locations using a 4th-order
462! Hamming time-stepping scheme.
463!-----------------------------------------------------------------------
464!
465 cff1=9.0_r8/8.0_r8
466 cff2=1.0_r8/8.0_r8
467 cff3=3.0_r8/8.0_r8
468 cff4=6.0_r8/8.0_r8
469 DO l=lstr,lend
470 IF (my_thread(l).and.bounded(l)) THEN
471 track(ixgrd,nfp1,l)=cff1*track(ixgrd,nf ,l)- &
472 & cff2*track(ixgrd,nfm2,l)+ &
473 & dt(ng)*(cff3*track(ixrhs,nfp1,l)+ &
474 & cff4*track(ixrhs,nf ,l)- &
475 & cff3*track(ixrhs,nfm1,l))
476 track(iygrd,nfp1,l)=cff1*track(iygrd,nf ,l)- &
477 & cff2*track(iygrd,nfm2,l)+ &
478 & dt(ng)*(cff3*track(iyrhs,nfp1,l)+ &
479 & cff4*track(iyrhs,nf ,l)- &
480 & cff3*track(iyrhs,nfm1,l))
481
482# if defined SOLVE3D && !defined FLOAT_VWALK
483!
484! Compute vertical position (grid units) 3D Lagrangian floats.
485!
486 IF (ftype(l).eq.flt_lagran) THEN
487# if defined FLOAT_BIOLOGY
488 track(izgrd,nfp1,l)=cff1*track(izgrd,nf ,l)- &
489 & cff2*track(izgrd,nfm2,l)+ &
490 & dt(ng)*(cff3*track(izrhs,nfp1,l)+ &
491 & cff4*track(izrhs,nf ,l)- &
492 & cff3*track(izrhs,nfm1,l)+ &
493 & cff3*track(iwbio,nfp1,l)* &
494 & track(i1ohz,nfp1,l)+ &
495 & cff4*track(iwbio,nf ,l)* &
496 & track(i1ohz,nf ,l)- &
497 & cff3*track(iwbio,nfm1,l)* &
498 & track(i1ohz,nf ,l))
499# else
500 track(izgrd,nfp1,l)=cff1*track(izgrd,nf ,l)- &
501 & cff2*track(izgrd,nfm2,l)+ &
502 & dt(ng)*(cff3*track(izrhs,nfp1,l)+ &
503 & cff4*track(izrhs,nf ,l)- &
504 & cff3*track(izrhs,nfm1,l))
505# endif
506!
507! Compute vertical position (grid units) for isobaric floats
508! (p=g*(z+zeta)=constant) or geopotential floats (constant depth).
509! Use bilinear interpolation to determine vertical position.
510!
511 ELSE IF ((ftype(l).eq.flt_isobar).or. &
512 & (ftype(l).eq.flt_geopot)) THEN
513 ir=int(track(ixgrd,nfp1,l))
514 jr=int(track(iygrd,nfp1,l))
515!
516 i1=min(max(ir ,0),lm(ng)+1)
517 i2=min(max(ir+1,1),lm(ng)+1)
518 j1=min(max(jr ,0),mm(ng)+1)
519 j2=min(max(jr+1,0),mm(ng)+1)
520!
521 p2=real(i2-i1,r8)*(track(ixgrd,nfp1,l)-real(i1,r8))
522 q2=real(j2-j1,r8)*(track(iygrd,nfp1,l)-real(j1,r8))
523 p1=1.0_r8-p2
524 q1=1.0_r8-q2
525# ifdef MASKING
526 cff7=p1*q1*grid(ng)%z_w(i1,j1,n(ng))*grid(ng)%rmask(i1,j1)+ &
527 & p2*q1*grid(ng)%z_w(i2,j1,n(ng))*grid(ng)%rmask(i2,j1)+ &
528 & p1*q2*grid(ng)%z_w(i1,j2,n(ng))*grid(ng)%rmask(i1,j2)+ &
529 & p2*q2*grid(ng)%z_w(i2,j2,n(ng))*grid(ng)%rmask(i2,j2)
530 cff8=p1*q1*grid(ng)%rmask(i1,j1)+ &
531 & p2*q1*grid(ng)%rmask(i2,j1)+ &
532 & p1*q2*grid(ng)%rmask(i1,j2)+ &
533 & p2*q2*grid(ng)%rmask(i2,j2)
534 IF (cff8.gt.0.0_r8) THEN
535 cff9=cff7/cff8
536 ELSE
537 cff9=0.0_r8
538 END IF
539# else
540 cff9=p1*q1*grid(ng)%z_w(i1,j1,n(ng))+ &
541 & p2*q1*grid(ng)%z_w(i2,j1,n(ng))+ &
542 & p1*q2*grid(ng)%z_w(i1,j2,n(ng))+ &
543 & p2*q2*grid(ng)%z_w(i2,j2,n(ng))
544# endif
545 cff6=cff9
546!
547 IF (ftype(l).eq.flt_geopot) THEN
548 zfloat=fz0(l)
549 ELSE IF (ftype(l).eq.flt_isobar) THEN
550 zfloat=fz0(l)+cff9
551 END IF
552!
553 DO k=n(ng)-1,0,-1
554# ifdef MASKING
555 cff7=p1*q1*grid(ng)%z_w(i1,j1,k)*grid(ng)%rmask(i1,j1)+ &
556 & p2*q1*grid(ng)%z_w(i2,j1,k)*grid(ng)%rmask(i2,j1)+ &
557 & p1*q2*grid(ng)%z_w(i1,j2,k)*grid(ng)%rmask(i1,j2)+ &
558 & p2*q2*grid(ng)%z_w(i2,j2,k)*grid(ng)%rmask(i2,j2)
559 cff8=p1*q1*grid(ng)%rmask(i1,j1)+ &
560 & p2*q1*grid(ng)%rmask(i2,j1)+ &
561 & p1*q2*grid(ng)%rmask(i1,j2)+ &
562 & p2*q2*grid(ng)%rmask(i2,j2)
563 cff5=0.0_r8
564 IF (cff8.gt.0.0_r8) cff5=cff7/cff8
565# else
566 cff5=p1*q1*grid(ng)%z_w(i1,j1,k)+ &
567 & p2*q1*grid(ng)%z_w(i2,j1,k)+ &
568 & p1*q2*grid(ng)%z_w(i1,j2,k)+ &
569 & p2*q2*grid(ng)%z_w(i2,j2,k)
570# endif
571 IF ((zfloat-cff5)*(cff6-zfloat).ge.0.0_r8) THEN
572 track(izgrd,nfp1,l)=real(k,r8)+(zfloat-cff5)/(cff6-cff5)
573 END IF
574 cff6=cff5
575 END DO
576 END IF
577# endif
578 END IF
579 END DO
580!
581!-----------------------------------------------------------------------
582! Determine floats status.
583!-----------------------------------------------------------------------
584!
585 IF (ewperiodic(ng)) THEN
586 cff1=real(lm(ng),r8)
587 DO l=lstr,lend
588 IF (my_thread(l).and.bounded(l)) THEN
589 IF (track(ixgrd,nfp1,l).ge.real(lm(ng)+1,r8)-0.5_r8) THEN
590 track(ixgrd,nfp1,l)=track(ixgrd,nfp1,l)-cff1
591 track(ixgrd,nf ,l)=track(ixgrd,nf ,l)-cff1
592 track(ixgrd,nfm1,l)=track(ixgrd,nfm1,l)-cff1
593 track(ixgrd,nfm2,l)=track(ixgrd,nfm2,l)-cff1
594 track(ixgrd,nfm3,l)=track(ixgrd,nfm3,l)-cff1
595 ELSE IF (track(ixgrd,nfp1,l).lt.0.5_r8) THEN
596 track(ixgrd,nfp1,l)=cff1+track(ixgrd,nfp1,l)
597 track(ixgrd,nf ,l)=cff1+track(ixgrd,nf ,l)
598 track(ixgrd,nfm1,l)=cff1+track(ixgrd,nfm1,l)
599 track(ixgrd,nfm2,l)=cff1+track(ixgrd,nfm2,l)
600 track(ixgrd,nfm3,l)=cff1+track(ixgrd,nfm3,l)
601 END IF
602 END IF
603 END DO
604# ifdef DISTRIBUTE
605 IF (ntilei(ng).gt.1) THEN
606 fwrk=reshape(track,(/npts/))
607 CALL mp_collect (ng, inlm, npts, fspv, fwrk)
608 track=reshape(fwrk,(/nfv(ng),nft+1,nfloats(ng)/))
609 DO l=lstr,lend
610 IF ((xstr.le.track(ixgrd,nfp1,l)).and. &
611 & (track(ixgrd,nfp1,l).lt.xend).and. &
612 & (ystr.le.track(iygrd,nfp1,l)).and. &
613 & (track(iygrd,nfp1,l).lt.yend)) THEN
614 my_thread(l)=.true.
615 ELSE IF (master.and.(.not.bounded(l))) THEN
616 my_thread(l)=.true.
617 ELSE
618 my_thread(l)=.false.
619 DO j=0,nft
620 DO i=1,nfv(ng)
621 track(i,j,l)=fspv
622 END DO
623 END DO
624 END IF
625 END DO
626 END IF
627# endif
628 ELSE
629 DO l=lstr,lend
630 IF (my_thread(l).and.bounded(l)) THEN
631 IF ((track(ixgrd,nfp1,l).ge.real(lm(ng)+1,r8)-0.5_r8).or. &
632 & (track(ixgrd,nfp1,l).lt.0.5_r8)) THEN
633 bounded(l)=.false.
634 END IF
635 END IF
636 END DO
637 END IF
638!
639 IF (nsperiodic(ng)) THEN
640 cff1=real(mm(ng),r8)
641 DO l=lstr,lend
642 IF (my_thread(l).and.bounded(l)) THEN
643 IF (track(iygrd,nfp1,l).ge.real(mm(ng)+1,r8)-0.5_r8) THEN
644 track(iygrd,nfp1,l)=track(iygrd,nfp1,l)-cff1
645 track(iygrd,nf ,l)=track(iygrd,nf ,l)-cff1
646 track(iygrd,nfm1,l)=track(iygrd,nfm1,l)-cff1
647 track(iygrd,nfm2,l)=track(iygrd,nfm2,l)-cff1
648 track(iygrd,nfm3,l)=track(iygrd,nfm3,l)-cff1
649 ELSE IF (track(iygrd,nfp1,l).lt.0.5_r8) THEN
650 track(iygrd,nfp1,l)=cff1+track(iygrd,nfp1,l)
651 track(iygrd,nf ,l)=cff1+track(iygrd,nf ,l)
652 track(iygrd,nfm1,l)=cff1+track(iygrd,nfm1,l)
653 track(iygrd,nfm2,l)=cff1+track(iygrd,nfm2,l)
654 track(iygrd,nfm3,l)=cff1+track(iygrd,nfm3,l)
655 END IF
656 END IF
657 END DO
658# ifdef DISTRIBUTE
659 IF (ntilej(ng).gt.1) THEN
660 fwrk=reshape(track,(/npts/))
661 CALL mp_collect (ng, inlm, npts, fspv, fwrk)
662 track=reshape(fwrk,(/nfv(ng),nft+1,nfloats(ng)/))
663 DO l=lstr,lend
664 IF ((xstr.le.track(ixgrd,nfp1,l)).and. &
665 & (track(ixgrd,nfp1,l).lt.xend).and. &
666 & (ystr.le.track(iygrd,nfp1,l)).and. &
667 & (track(iygrd,nfp1,l).lt.yend)) THEN
668 my_thread(l)=.true.
669 ELSE IF (master.and.(.not.bounded(l))) THEN
670 my_thread(l)=.true.
671 ELSE
672 my_thread(l)=.false.
673 DO j=0,nft
674 DO i=1,nfv(ng)
675 track(i,j,l)=fspv
676 END DO
677 END DO
678 END IF
679 END DO
680 END IF
681# endif
682 ELSE
683 DO l=lstr,lend
684 IF (my_thread(l).and.bounded(l)) THEN
685 IF ((track(iygrd,nfp1,l).ge.real(mm(ng)+1,r8)-0.5_r8).or. &
686 & (track(iygrd,nfp1,l).lt.0.5_r8)) THEN
687 bounded(l)=.false.
688 END IF
689 END IF
690 END DO
691 END IF
692!
693!-----------------------------------------------------------------------
694! If appropriate, activate the release of new floats and set initial
695! positions for all time levels.
696!-----------------------------------------------------------------------
697!
698 halfdt=0.5_r8*dt(ng)
699
700 DO l=lstr,lend
701 IF (.not.bounded(l).and. &
702 & (time(ng)-halfdt.le.tinfo(itstr,l).and. &
703 & time(ng)+halfdt.gt.tinfo(itstr,l))) THEN
704 bounded(l)=.true.
705 IF ((tinfo(ixgrd,l).lt.0.5_r8).or. &
706 & (tinfo(iygrd,l).lt.0.5_r8).or. &
707 & (tinfo(ixgrd,l).gt.real(lm(ng),r8)+0.5_r8).or. &
708 & (tinfo(iygrd,l).gt.real(mm(ng),r8)+0.5_r8)) THEN
709 bounded(l)=.false. ! outside application grid
710 END IF
711# if defined SOLVE3D && defined FLOAT_STICKY
712 stuck(l)=.false.
713# endif
714# ifdef DISTRIBUTE
715 IF ((xstr.le.tinfo(ixgrd,l)).and. &
716 & (tinfo(ixgrd,l).lt.xend).and. &
717 & (ystr.le.tinfo(iygrd,l)).and. &
718 & (tinfo(iygrd,l).lt.yend).and.bounded(l)) THEN
719 DO j=0,nft
720 track(ixgrd,j,l)=tinfo(ixgrd,l)
721 track(iygrd,j,l)=tinfo(iygrd,l)
722# ifdef SOLVE3D
723 track(izgrd,j,l)=tinfo(izgrd,l)
724# endif
725 END DO
726 my_thread(l)=.true.
727 ELSE
728 my_thread(l)=.false.
729 DO j=0,nft
730 DO i=1,nfv(ng)
731 track(i,j,l)=fspv
732 END DO
733 END DO
734 END IF
735# else
736 IF (bounded(l)) THEN
737 my_thread(l)=.true.
738 DO j=0,nft
739 track(ixgrd,j,l)=tinfo(ixgrd,l)
740 track(iygrd,j,l)=tinfo(iygrd,l)
741# ifdef SOLVE3D
742 track(izgrd,j,l)=tinfo(izgrd,l)
743# endif
744 END DO
745 END IF
746# endif
747 END IF
748 END DO
749!
750!-----------------------------------------------------------------------
751! Calculate slopes with corrected locations.
752!-----------------------------------------------------------------------
753!
754# ifdef SOLVE3D
755 CALL interp_floats (ng, lbi, ubi, lbj, ubj, 1, n(ng), &
756 & lstr, lend, nfp1, ixrhs, isuvel, &
757 & -u3dvar, lmask, spval, nudg, &
758 & grid(ng) % pm, &
759 & grid(ng) % pn, &
760 & grid(ng) % Hz, &
761# ifdef MASKING
762 & grid(ng) % rmask, &
763# endif
764 & ocean(ng) % u(:,:,:,nnew), &
765 & my_thread, bounded, track)
766
767 CALL interp_floats (ng, lbi, ubi, lbj, ubj, 1, n(ng), &
768 & lstr, lend, nfp1, iyrhs, isvvel, &
769 & -v3dvar, lmask, spval, nudg, &
770 & grid(ng) % pm, &
771 & grid(ng) % pn, &
772 & grid(ng) % Hz, &
773# ifdef MASKING
774 & grid(ng) % rmask, &
775# endif
776 & ocean(ng) % v(:,:,:,nnew), &
777 & my_thread, bounded, track)
778
779# if !defined FLOAT_VWALK
780 CALL interp_floats (ng, lbi, ubi, lbj, ubj, 0, n(ng), &
781 & lstr, lend, nfp1, izrhs, isbw3d, &
782 & -w3dvar, lmask, spval, nudg, &
783 & grid(ng) % pm, &
784 & grid(ng) % pn, &
785 & grid(ng) % Hz, &
786# ifdef MASKING
787 & grid(ng) % rmask, &
788# endif
789 & ocean(ng) % W, &
790 & my_thread, bounded, track)
791# endif
792# ifdef FLOAT_BIOLOGY
793 CALL interp_floats (ng, lbi, ubi, lbj, ubj, 1, n(ng), &
794 & lstr, lend, nfp1, i1ohz, isbr3d, &
795 & r3dvar, gmask, fspv, nudg, &
796 & grid(ng) % pm, &
797 & grid(ng) % pn, &
798 & grid(ng) % Hz, &
799# ifdef MASKING
800 & grid(ng) % rmask, &
801# endif
802 & grid(ng) % Hz, &
803 & my_thread, bounded, track)
804 DO l=lstr,lend
805 IF (my_thread(l).and.bounded(l)) THEN
806 track(i1ohz,nfp1,l)=1.0_r8/track(i1ohz,nfp1,l)
807 END IF
808 END DO
809# endif
810# else
811 CALL interp_floats (ng, lbi, ubi, lbj, ubj, 1, 1, &
812 & lstr, lend, nfp1, ixrhs, isubar, &
813 & -u2dvar, lmask, spval, nudg, &
814 & grid(ng) % pm, &
815 & grid(ng) % pn, &
816# ifdef MASKING
817 & grid(ng) % rmask, &
818# endif
819 & ocean(ng) % ubar(:,:,knew), &
820 & my_thread, bounded, track)
821
822 CALL interp_floats (ng, lbi, ubi, lbj, ubj, 1, 1, &
823 & lstr, lend, nfp1, iyrhs, isvbar, &
824 & -v2dvar, lmask, spval, nudg, &
825 & grid(ng) % pm, &
826 & grid(ng) % pn, &
827# ifdef MASKING
828 & grid(ng) % rmask, &
829# endif
830 & ocean(ng) % vbar(:,:,knew), &
831 & my_thread, bounded, track)
832# endif
833!
834! If newly released floats, initialize slopes at all time levels.
835!
836 DO l=lstr,lend
837 IF (my_thread(l).and.bounded(l).and. &
838 & (time(ng)-halfdt.le.tinfo(itstr,l).and. &
839 & time(ng)+halfdt.gt.tinfo(itstr,l))) THEN
840 xrhs=track(ixrhs,nfp1,l)
841 yrhs=track(iyrhs,nfp1,l)
842# ifdef SOLVE3D
843 zrhs=track(izrhs,nfp1,l)
844# endif
845# ifdef FLOAT_BIOLOGY
846 ohz=track(i1ohz,nfp1,l)
847# endif
848 DO i=0,nft
849 track(ixrhs,i,l)=xrhs
850 track(iyrhs,i,l)=yrhs
851# ifdef SOLVE3D
852 track(izrhs,i,l)=zrhs
853# endif
854# ifdef FLOAT_BIOLOGY
855 track(i1ohz,i,l)=ohz
856# endif
857 END DO
858 END IF
859 END DO
860!
861!-----------------------------------------------------------------------
862! Interpolate various output variables at the corrected locations.
863!-----------------------------------------------------------------------
864!
865 IF (spherical) THEN
866 CALL interp_floats (ng, lbi, ubi, lbj, ubj, 1, 1, &
867 & lstr, lend, nfp1, iflon, isbr2d, &
868 & r2dvar, gmask, spval, nudg, &
869 & grid(ng) % pm, &
870 & grid(ng) % pn, &
871# ifdef SOLVE3D
872 & grid(ng) % Hz, &
873# endif
874# ifdef MASKING
875 & grid(ng) % rmask, &
876# endif
877 & grid(ng) % lonr, &
878 & my_thread, bounded, track)
879
880 CALL interp_floats (ng, lbi, ubi, lbj, ubj, 1, 1, &
881 & lstr, lend, nfp1, iflat, isbr2d, &
882 & r2dvar, gmask, spval, nudg, &
883 & grid(ng) % pm, &
884 & grid(ng) % pn, &
885# ifdef SOLVE3D
886 & grid(ng) % Hz, &
887# endif
888# ifdef MASKING
889 & grid(ng) % rmask, &
890# endif
891 & grid(ng) % latr, &
892 & my_thread, bounded, track)
893 ELSE
894 CALL interp_floats (ng, lbi, ubi, lbj, ubj, 1, 1, &
895 & lstr, lend, nfp1, iflon, isbr2d, &
896 & r2dvar, gmask, spval, nudg, &
897 & grid(ng) % pm, &
898 & grid(ng) % pn, &
899# ifdef SOLVE3D
900 & grid(ng) % Hz, &
901# endif
902# ifdef MASKING
903 & grid(ng) % rmask, &
904# endif
905 & grid(ng) % xr, &
906 & my_thread, bounded, track)
907
908 CALL interp_floats (ng, lbi, ubi, lbj, ubj, 1, 1, &
909 & lstr, lend, nfp1, iflat, isbr2d, &
910 & r2dvar, gmask, spval, nudg, &
911 & grid(ng) % pm, &
912 & grid(ng) % pn, &
913# ifdef SOLVE3D
914 & grid(ng) % Hz, &
915# endif
916# ifdef MASKING
917 & grid(ng) % rmask, &
918# endif
919 & grid(ng) % yr, &
920 & my_thread, bounded, track)
921 END IF
922# ifdef SOLVE3D
923 CALL interp_floats (ng, lbi, ubi, lbj, ubj, 0, n(ng), &
924 & lstr, lend, nfp1, idpth, isbw3d, &
925 & w3dvar, lmask, spval, nudg, &
926 & grid(ng) % pm, &
927 & grid(ng) % pn, &
928 & grid(ng) % Hz, &
929# ifdef MASKING
930 & grid(ng) % rmask, &
931# endif
932 & grid(ng) % z_w, &
933 & my_thread, bounded, track)
934
935 CALL interp_floats (ng, lbi, ubi, lbj, ubj, 1, n(ng), &
936 & lstr, lend, nfp1, ifden, isbr3d, &
937 & r3dvar, lmask, spval, nudg, &
938 & grid(ng) % pm, &
939 & grid(ng) % pn, &
940 & grid(ng) % Hz, &
941# ifdef MASKING
942 & grid(ng) % rmask, &
943# endif
944 & ocean(ng) % rho, &
945 & my_thread, bounded, track)
946
947 DO itrc=1,nt(ng)
948 CALL interp_floats (ng, lbi, ubi, lbj, ubj, 1, n(ng), &
949 & lstr, lend, nfp1, iftvar(itrc), &
950 & istvar(itrc), r3dvar, lmask, spval, nudg, &
951 & grid(ng) % pm, &
952 & grid(ng) % pn, &
953 & grid(ng) % Hz, &
954# ifdef MASKING
955 & grid(ng) % rmask, &
956# endif
957 & ocean(ng) % t(:,:,:,nnew,itrc), &
958 & my_thread, bounded, track)
959 END DO
960# endif
961# ifdef FLOAT_BIOLOGY
962!
963!-----------------------------------------------------------------------
964! Compute biological float behavior, corrector step.
965!-----------------------------------------------------------------------
966!
967 CALL biology_floats (ng, lstr, lend, .false., my_thread)
968# endif
969# if defined SOLVE3D && defined FLOAT_VWALK && !defined VWALK_FORWARD
970!
971!-----------------------------------------------------------------------
972! Compute vertical positions due to vertical random walk, corrector
973! step.
974!-----------------------------------------------------------------------
975!
976 CALL vwalk_floats (ng, lstr, lend, .false., my_thread, nudg)
977# endif
978# ifdef SOLVE3D
979# ifdef FLOAT_STICKY
980!
981! Floats that hit the surface are reflected; floats that hit the bottom
982! get stuck.
983!
984 DO l=lstr,lend
985 IF (my_thread(l).and.bounded(l)) THEN
986 IF (stuck(l)) THEN
987 track(ixgrd,nfp1,l)=track(ixgrd,nf,l)
988 track(iygrd,nfp1,l)=track(iygrd,nf,l)
989 track(izgrd,nfp1,l)=track(izgrd,nf,l)
990 ELSE
991 IF (track(izgrd,nfp1,l).gt.real(n(ng),r8)) THEN
992 DO j=0,nft
993 track(izgrd,j,l)=2.0_r8*real(n(ng),r8)- &
994 & track(izgrd,j,l)
995 END DO
996 ELSE IF (track(izgrd,nfp1,l).lt.0.0_r8) THEN
997 DO j=0,nft
998 track(izgrd,j,l)=0.0_r8
999 END DO
1000 stuck(l)=.true.
1001 END IF
1002 END IF
1003 END IF
1004 END DO
1005# else
1006!
1007! Floats that hit the surface or bottom are reflected
1008!
1009 DO l=lstr,lend
1010 IF (my_thread(l).and.bounded(l)) THEN
1011 IF (track(izgrd,nfp1,l).gt.real(n(ng),r8)) THEN
1012 DO j=0,nft
1013 track(izgrd,j,l)=2.0_r8*real(n(ng),r8)-track(izgrd,j,l)
1014 END DO
1015 ELSE IF (track(izgrd,nfp1,l).lt.0.0_r8) THEN
1016 DO j=0,nft
1017 track(izgrd,:,l)=-track(izgrd,:,l)
1018 END DO
1019 END IF
1020 END IF
1021 END DO
1022# endif
1023# endif
1024# ifdef DISTRIBUTE
1025!
1026!-----------------------------------------------------------------------
1027! Collect floats on all nodes.
1028!-----------------------------------------------------------------------
1029!
1030 fwrk=reshape(track,(/npts/))
1031 CALL mp_collect (ng, inlm, npts, fspv, fwrk)
1032 track=reshape(fwrk,(/nfv(ng),nft+1,nfloats(ng)/))
1033!
1034! Collect the bounded status switch.
1035!
1036 fwrk=fspv
1037 DO l=1,nfloats(ng)
1038 IF (bounded(l)) THEN
1039 fwrk(l)=1.0_r8
1040 END IF
1041 END DO
1042 CALL mp_collect (ng, inlm, nfloats(ng), fspv, fwrk)
1043 DO l=1,nfloats(ng)
1044 IF (fwrk(l).ne.fspv) THEN
1045 bounded(l)=.true.
1046 ELSE
1047 bounded(l)=.false.
1048 END IF
1049 END DO
1050# endif
1051!
1052 RETURN
subroutine, public biology_floats(ng, lstr, lend, predictor, my_thread)
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 flt_isobar
Definition mod_floats.F:126
integer, parameter iflat
Definition mod_floats.F:85
integer, parameter i1ohz
Definition mod_floats.F:97
integer, parameter ifden
Definition mod_floats.F:90
integer, parameter iflon
Definition mod_floats.F:84
integer, parameter idpth
Definition mod_floats.F:86
integer, parameter iwbio
Definition mod_floats.F:101
integer, parameter flt_geopot
Definition mod_floats.F:127
integer, parameter iyrhs
Definition mod_floats.F:88
integer, parameter iygrd
Definition mod_floats.F:82
integer, dimension(:), allocatable iftvar
Definition mod_floats.F:116
integer, parameter itstr
Definition mod_floats.F:80
integer, parameter izrhs
Definition mod_floats.F:89
integer, parameter ixgrd
Definition mod_floats.F:81
integer, parameter izgrd
Definition mod_floats.F:83
integer, parameter ixrhs
Definition mod_floats.F:87
integer, parameter flt_lagran
Definition mod_floats.F:125
type(t_grid), dimension(:), allocatable grid
Definition mod_grid.F:365
integer isvvel
integer isvbar
integer, dimension(:), allocatable istvar
integer isuvel
integer isubar
integer isbw3d
integer isbr3d
integer isbr2d
type(t_ocean), dimension(:), allocatable ocean
Definition mod_ocean.F:351
logical master
integer, dimension(:), allocatable nfv
Definition mod_param.F:547
integer, dimension(:), allocatable nfloats
Definition mod_param.F:543
integer, dimension(:), allocatable n
Definition mod_param.F:479
type(t_bounds), dimension(:), allocatable bounds
Definition mod_param.F:232
integer, parameter r3dvar
Definition mod_param.F:721
integer, parameter u3dvar
Definition mod_param.F:722
integer, parameter nft
Definition mod_param.F:539
integer, dimension(:), allocatable lm
Definition mod_param.F:455
integer, parameter u2dvar
Definition mod_param.F:718
integer, dimension(:), allocatable ntilei
Definition mod_param.F:677
integer, parameter w3dvar
Definition mod_param.F:724
integer, dimension(:), allocatable nt
Definition mod_param.F:489
integer, dimension(:), allocatable mm
Definition mod_param.F:456
integer, parameter r2dvar
Definition mod_param.F:717
integer, parameter v2dvar
Definition mod_param.F:719
integer, parameter v3dvar
Definition mod_param.F:723
integer, dimension(:), allocatable ntilej
Definition mod_param.F:678
logical spherical
real(dp), parameter spval
real(dp), dimension(:), allocatable dt
logical, dimension(:), allocatable ewperiodic
logical, dimension(:), allocatable nsperiodic
real(dp), dimension(:), allocatable time
subroutine, public vwalk_floats(ng, lstr, lend, predictor, my_thread, nudg)

References biology_floats_mod::biology_floats(), mod_param::bounds, mod_scalars::dt, mod_scalars::ewperiodic, mod_floats::flt_geopot, mod_floats::flt_isobar, mod_floats::flt_lagran, mod_grid::grid, mod_floats::i1ohz, mod_floats::idpth, mod_floats::ifden, mod_floats::iflat, mod_floats::iflon, mod_floats::iftvar, mod_param::inlm, interp_floats_mod::interp_floats(), mod_ncparam::isbr2d, mod_ncparam::isbr3d, mod_ncparam::isbw3d, mod_ncparam::istvar, mod_ncparam::isubar, mod_ncparam::isuvel, mod_ncparam::isvbar, mod_ncparam::isvvel, mod_floats::itstr, mod_floats::iwbio, mod_floats::ixgrd, mod_floats::ixrhs, mod_floats::iygrd, mod_floats::iyrhs, mod_floats::izgrd, mod_param::lm, mod_parallel::master, mod_param::mm, mod_parallel::myrank, mod_param::n, mod_scalars::nsperiodic, mod_param::nt, mod_param::ntilei, mod_param::ntilej, mod_ocean::ocean, mod_param::r2dvar, mod_param::r3dvar, mod_scalars::spherical, mod_scalars::spval, mod_scalars::time, mod_param::u2dvar, mod_param::u3dvar, mod_param::v2dvar, mod_param::v3dvar, vwalk_floats_mod::vwalk_floats(), and mod_param::w3dvar.

Referenced by step_floats().

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