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

Functions/Subroutines

subroutine, public wrt_floats (ng)
 
subroutine, private wrt_floats_nf90 (ng)
 
subroutine, private wrt_floats_pio (ng)
 

Function/Subroutine Documentation

◆ wrt_floats()

subroutine, public wrt_floats_mod::wrt_floats ( integer, intent(in) ng)

Definition at line 39 of file wrt_floats.F.

40!***********************************************************************
41!
42! Imported variable declarations.
43!
44 integer, intent(in) :: ng
45!
46! Local variable declarations.
47!
48 character (len=*), parameter :: MyFile = &
49 & __FILE__
50!
51!-----------------------------------------------------------------------
52! Write out history fields according to IO type.
53!-----------------------------------------------------------------------
54!
55 SELECT CASE (flt(ng)%IOtype)
56 CASE (io_nf90)
57 CALL wrt_floats_nf90 (ng)
58
59# if defined PIO_LIB && defined DISTRIBUTE
60 CASE (io_pio)
61 CALL wrt_floats_pio (ng)
62# endif
63 CASE DEFAULT
64 IF (master) WRITE (stdout,10) flt(ng)%IOtype
65 exit_flag=3
66 END SELECT
67 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
68!
69 10 FORMAT (' WRT_FLOATS - Illegal output file type, io_type = ',i0, &
70 & /,14x,'Check KeyWord ''OUT_LIB'' in ''roms.in''.')
71!
72 RETURN

References mod_scalars::exit_flag, mod_iounits::flt, strings_mod::founderror(), mod_ncparam::io_nf90, mod_ncparam::io_pio, mod_parallel::master, mod_scalars::noerror, mod_iounits::stdout, wrt_floats_nf90(), and wrt_floats_pio().

Referenced by output().

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

◆ wrt_floats_nf90()

subroutine, private wrt_floats_mod::wrt_floats_nf90 ( integer, intent(in) ng)
private

Definition at line 76 of file wrt_floats.F.

77!***********************************************************************
78!
79 USE mod_netcdf
80!
81! Imported variable declarations.
82!
83 integer, intent(in) :: ng
84!
85! Local variable declarations.
86!
87 integer :: Fcount, itrc, l, status
88
89 real(r8), dimension(Nfloats(ng)) :: Tout
90!
91 character (len=*), parameter :: MyFile = &
92 & __FILE__//", wrt_floats_nf90"
93!
94 sourcefile=myfile
95!
96!-----------------------------------------------------------------------
97! Write out station data at RHO-points.
98!-----------------------------------------------------------------------
99!
100 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
101!
102! Set time record index.
103!
104 flt(ng)%Rindex=flt(ng)%Rindex+1
105 fcount=flt(ng)%Fcount
106 flt(ng)%Nrec(fcount)=flt(ng)%Nrec(fcount)+1
107!
108! Write out model time (s).
109!
110 CALL netcdf_put_fvar (ng, inlm, flt(ng)%name, &
111 & trim(vname(1,idtime)), time(ng:), &
112 & (/flt(ng)%Rindex/), (/1/), &
113 & ncid = flt(ng)%ncid, &
114 & varid = flt(ng)%Vid(idtime))
115 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
116!
117! Write out floats X-grid locations.
118!
119 DO l=1,nfloats(ng)
120 IF (drifter(ng)%bounded(l)) THEN
121 tout(l)=drifter(ng)%track(ixgrd,nf(ng),l)
122 ELSE
123 tout(l)=spval
124 END IF
125 END DO
126 CALL netcdf_put_fvar (ng, inlm, flt(ng)%name, &
127 & 'Xgrid', tout, &
128 & (/1,flt(ng)%Rindex/), (/nfloats(ng),1/), &
129 & ncid = flt(ng)%ncid, &
130 & varid = flt(ng)%Vid(idxgrd))
131 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
132!
133! Write out floats Y-grid locations.
134!
135 DO l=1,nfloats(ng)
136 IF (drifter(ng)%bounded(l)) THEN
137 tout(l)=drifter(ng)%track(iygrd,nf(ng),l)
138 ELSE
139 tout(l)=spval
140 END IF
141 END DO
142 CALL netcdf_put_fvar (ng, inlm, flt(ng)%name, &
143 & 'Ygrid', tout, &
144 & (/1,flt(ng)%Rindex/), (/nfloats(ng),1/), &
145 & ncid = flt(ng)%ncid, &
146 & varid = flt(ng)%Vid(idygrd))
147 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
148
149# ifdef SOLVE3D
150!
151! Write out floats Z-grid locations.
152!
153 DO l=1,nfloats(ng)
154 IF (drifter(ng)%bounded(l)) THEN
155 tout(l)=drifter(ng)%track(izgrd,nf(ng),l)
156 ELSE
157 tout(l)=spval
158 END IF
159 END DO
160 CALL netcdf_put_fvar (ng, inlm, flt(ng)%name, &
161 & 'Zgrid', tout, &
162 & (/1,flt(ng)%Rindex/), (/nfloats(ng),1/), &
163 & ncid = flt(ng)%ncid, &
164 & varid = flt(ng)%Vid(idzgrd))
165 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
166# endif
167!
168! Write out floats (lon,lat) or (x,y) locations.
169!
170 DO l=1,nfloats(ng)
171 IF (drifter(ng)%bounded(l)) THEN
172 tout(l)=drifter(ng)%track(iflon,nf(ng),l)
173 ELSE
174 tout(l)=spval
175 END IF
176 END DO
177 IF (spherical) THEN
178 CALL netcdf_put_fvar (ng, inlm, flt(ng)%name, &
179 & 'lon', tout, &
180 & (/1,flt(ng)%Rindex/), (/nfloats(ng),1/), &
181 & ncid = flt(ng)%ncid, &
182 & varid = flt(ng)%Vid(idglon))
183 ELSE
184 CALL netcdf_put_fvar (ng, inlm, flt(ng)%name, &
185 & 'x', tout, &
186 & (/1,flt(ng)%Rindex/), (/nfloats(ng),1/), &
187 & ncid = flt(ng)%ncid, &
188 & varid = flt(ng)%Vid(idglon))
189 END IF
190 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
191!
192 DO l=1,nfloats(ng)
193 IF (drifter(ng)%bounded(l)) THEN
194 tout(l)=drifter(ng)%track(iflat,nf(ng),l)
195 ELSE
196 tout(l)=spval
197 END IF
198 END DO
199 IF (spherical) THEN
200 CALL netcdf_put_fvar (ng, inlm, flt(ng)%name, &
201 & 'lat', tout, &
202 & (/1,flt(ng)%Rindex/), (/nfloats(ng),1/), &
203 & ncid = flt(ng)%ncid, &
204 & varid = flt(ng)%Vid(idglat))
205 ELSE
206 CALL netcdf_put_fvar (ng, inlm, flt(ng)%name, &
207 & 'y', tout, &
208 & (/1,flt(ng)%Rindex/), (/nfloats(ng),1/), &
209 & ncid = flt(ng)%ncid, &
210 & varid = flt(ng)%Vid(idglat))
211 END IF
212 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
213
214# ifdef SOLVE3D
215!
216! Write out floats depths.
217!
218 DO l=1,nfloats(ng)
219 IF (drifter(ng)%bounded(l)) THEN
220 tout(l)=drifter(ng)%track(idpth,nf(ng),l)
221 ELSE
222 tout(l)=spval
223 END IF
224 END DO
225 CALL netcdf_put_fvar (ng, inlm, flt(ng)%name, &
226 & 'depth', tout, &
227 & (/1,flt(ng)%Rindex/), (/nfloats(ng),1/), &
228 & ncid = flt(ng)%ncid, &
229 & varid = flt(ng)%Vid(iddpth))
230 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
231!
232! Write out density anomaly.
233!
234 DO l=1,nfloats(ng)
235 IF (drifter(ng)%bounded(l)) THEN
236 tout(l)=drifter(ng)%track(ifden,nf(ng),l)
237 ELSE
238 tout(l)=spval
239 END IF
240 END DO
241 CALL netcdf_put_fvar (ng, inlm, flt(ng)%name, &
242 & trim(vname(1,iddano)), tout, &
243 & (/1,flt(ng)%Rindex/), (/nfloats(ng),1/), &
244 & ncid = flt(ng)%ncid, &
245 & varid = flt(ng)%Vid(iddano))
246 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
247!
248! Write out tracer type variables.
249!
250 DO itrc=1,nt(ng)
251 DO l=1,nfloats(ng)
252 IF (drifter(ng)%bounded(l)) THEN
253 tout(l)=drifter(ng)%track(iftvar(itrc),nf(ng),l)
254 ELSE
255 tout(l)=spval
256 END IF
257 END DO
258 CALL netcdf_put_fvar (ng, inlm, flt(ng)%name, &
259 & trim(vname(1,idtvar(itrc))), tout, &
260 & (/1,flt(ng)%Rindex/), (/nfloats(ng),1/), &
261 & ncid = flt(ng)%ncid, &
262 & varid = flt(ng)%Tid(itrc))
263 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
264 END DO
265# endif
266
267# ifdef FLOAT_OYSTER
268!
269! Write out biological float swimming time.
270!
271 DO l=1,nfloats(ng)
272 IF (drifter(ng)%bounded(l)) THEN
273 tout(l)=drifter(ng)%track(iswim,nf(ng),l)
274 ELSE
275 tout(l)=spval
276 END IF
277 END DO
278 CALL netcdf_put_fvar (ng, inlm, flt(ng)%name, &
279 & 'swim_time', tout, &
280 & (/1,flt(ng)%Rindex/), (/nfloats(ng),1/), &
281 & ncid = flt(ng)%ncid, &
282 & varid = flt(ng)%Vid(idswim))
283 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
284!
285! Write out biological float vertical velocity.
286!
287 DO l=1,nfloats(ng)
288 IF (drifter(ng)%bounded(l)) THEN
289 tout(l)=drifter(ng)%track(iwbio,nf(ng),l)
290 ELSE
291 tout(l)=spval
292 END IF
293 END DO
294 CALL netcdf_put_fvar (ng, inlm, flt(ng)%name, &
295 & 'w_bio', tout, &
296 & (/1,flt(ng)%Rindex/), (/nfloats(ng),1/), &
297 & ncid = flt(ng)%ncid, &
298 & varid = flt(ng)%Vid(idwbio))
299 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
300!
301! Write out biological float size (length).
302!
303 DO l=1,nfloats(ng)
304 IF (drifter(ng)%bounded(l)) THEN
305 tout(l)=drifter(ng)%track(isizf,nf(ng),l)
306 ELSE
307 tout(l)=spval
308 END IF
309 END DO
310 CALL netcdf_put_fvar (ng, inlm, flt(ng)%name, &
311 & 'bio_size', tout, &
312 & (/1,flt(ng)%Rindex/), (/nfloats(ng),1/), &
313 & ncid = flt(ng)%ncid, &
314 & varid = flt(ng)%Vid(idsize))
315 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
316!
317! Write out biological float sinking velocity.
318!
319 DO l=1,nfloats(ng)
320 IF (drifter(ng)%bounded(l)) THEN
321 tout(l)=drifter(ng)%track(iwsin,nf(ng),l)
322 ELSE
323 tout(l)=spval
324 END IF
325 END DO
326 CALL netcdf_put_fvar (ng, inlm, flt(ng)%name, &
327 & 'bio_sink', tout, &
328 & (/1,flt(ng)%Rindex/), (/nfloats(ng),1/), &
329 & ncid = flt(ng)%ncid, &
330 & varid = flt(ng)%Vid(idwsin))
331 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
332# endif
333!
334!-----------------------------------------------------------------------
335! Synchronize floats NetCDF file to disk.
336!-----------------------------------------------------------------------
337!
338 CALL netcdf_sync (ng, inlm, flt(ng)%name, flt(ng)%ncid)
339!
340 RETURN
subroutine, public netcdf_sync(ng, model, ncname, ncid)

References mod_floats::drifter, mod_scalars::exit_flag, mod_iounits::flt, strings_mod::founderror(), mod_ncparam::iddano, mod_ncparam::iddpth, mod_ncparam::idglat, mod_ncparam::idglon, mod_floats::idpth, mod_ncparam::idsize, mod_ncparam::idswim, mod_ncparam::idtime, mod_ncparam::idtvar, mod_ncparam::idwbio, mod_ncparam::idwsin, mod_ncparam::idxgrd, mod_ncparam::idygrd, mod_ncparam::idzgrd, mod_floats::ifden, mod_floats::iflat, mod_floats::iflon, mod_floats::iftvar, mod_param::inlm, mod_floats::isizf, mod_floats::iswim, mod_floats::iwbio, mod_floats::iwsin, mod_floats::ixgrd, mod_floats::iygrd, mod_floats::izgrd, mod_netcdf::netcdf_sync(), mod_stepping::nf, mod_param::nfloats, mod_scalars::noerror, mod_param::nt, mod_iounits::sourcefile, mod_scalars::spherical, mod_scalars::spval, mod_scalars::time, and mod_ncparam::vname.

Referenced by wrt_floats().

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

◆ wrt_floats_pio()

subroutine, private wrt_floats_mod::wrt_floats_pio ( integer, intent(in) ng)
private

Definition at line 346 of file wrt_floats.F.

347!***********************************************************************
348!
350!
351! Imported variable declarations.
352!
353 integer, intent(in) :: ng
354!
355! Local variable declarations.
356!
357 integer :: Fcount, itrc, l, status
358
359 real(r8), dimension(Nfloats(ng)) :: Tout
360!
361 character (len=*), parameter :: MyFile = &
362 & __FILE__//", wrt_floats_pio"
363!
364 sourcefile=myfile
365!
366!-----------------------------------------------------------------------
367! Write out station data at RHO-points.
368!-----------------------------------------------------------------------
369!
370 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
371!
372! Set time record index.
373!
374 flt(ng)%Rindex=flt(ng)%Rindex+1
375 fcount=flt(ng)%Fcount
376 flt(ng)%Nrec(fcount)=flt(ng)%Nrec(fcount)+1
377!
378! Write out model time (s).
379!
380 CALL pio_netcdf_put_fvar (ng, inlm, flt(ng)%name, &
381 & trim(vname(1,idtime)), time(ng:), &
382 & (/flt(ng)%Rindex/), (/1/), &
383 & piofile = flt(ng)%pioFile, &
384 & piovar = flt(ng)%pioVar(idtime)%vd)
385 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
386!
387! Write out floats X-grid locations.
388!
389 DO l=1,nfloats(ng)
390 IF (drifter(ng)%bounded(l)) THEN
391 tout(l)=drifter(ng)%track(ixgrd,nf(ng),l)
392 ELSE
393 tout(l)=spval
394 END IF
395 END DO
396 CALL pio_netcdf_put_fvar (ng, inlm, flt(ng)%name, &
397 & 'Xgrid', tout, &
398 & (/1,flt(ng)%Rindex/), &
399 & (/nfloats(ng),1/), &
400 & piofile = flt(ng)%pioFile, &
401 & piovar = flt(ng)%pioVar(idxgrd)%vd)
402 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
403!
404! Write out floats Y-grid locations.
405!
406 DO l=1,nfloats(ng)
407 IF (drifter(ng)%bounded(l)) THEN
408 tout(l)=drifter(ng)%track(iygrd,nf(ng),l)
409 ELSE
410 tout(l)=spval
411 END IF
412 END DO
413 CALL pio_netcdf_put_fvar (ng, inlm, flt(ng)%name, &
414 & 'Ygrid', tout, &
415 & (/1,flt(ng)%Rindex/), &
416 & (/nfloats(ng),1/), &
417 & piofile = flt(ng)%pioFile, &
418 & piovar = flt(ng)%pioVar(idygrd)%vd)
419 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
420
421# ifdef SOLVE3D
422!
423! Write out floats Z-grid locations.
424!
425 DO l=1,nfloats(ng)
426 IF (drifter(ng)%bounded(l)) THEN
427 tout(l)=drifter(ng)%track(izgrd,nf(ng),l)
428 ELSE
429 tout(l)=spval
430 END IF
431 END DO
432 CALL pio_netcdf_put_fvar (ng, inlm, flt(ng)%name, &
433 & 'Zgrid', tout, &
434 & (/1,flt(ng)%Rindex/), &
435 & (/nfloats(ng),1/), &
436 & piofile = flt(ng)%pioFile, &
437 & piovar = flt(ng)%pioVar(idzgrd)%vd)
438 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
439# endif
440!
441! Write out floats (lon,lat) or (x,y) locations.
442!
443 DO l=1,nfloats(ng)
444 IF (drifter(ng)%bounded(l)) THEN
445 tout(l)=drifter(ng)%track(iflon,nf(ng),l)
446 ELSE
447 tout(l)=spval
448 END IF
449 END DO
450 IF (spherical) THEN
451 CALL pio_netcdf_put_fvar (ng, inlm, flt(ng)%name, &
452 & 'lon', tout, &
453 & (/1,flt(ng)%Rindex/), &
454 & (/nfloats(ng),1/), &
455 & piofile = flt(ng)%pioFile, &
456 & piovar = flt(ng)%pioVar(idglon)%vd)
457 ELSE
458 CALL pio_netcdf_put_fvar (ng, inlm, flt(ng)%name, &
459 & 'x', tout, &
460 & (/1,flt(ng)%Rindex/), &
461 & (/nfloats(ng),1/), &
462 & piofile = flt(ng)%pioFile, &
463 & piovar = flt(ng)%pioVar(idglon)%vd)
464 END IF
465 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
466!
467 DO l=1,nfloats(ng)
468 IF (drifter(ng)%bounded(l)) THEN
469 tout(l)=drifter(ng)%track(iflat,nf(ng),l)
470 ELSE
471 tout(l)=spval
472 END IF
473 END DO
474 IF (spherical) THEN
475 CALL pio_netcdf_put_fvar (ng, inlm, flt(ng)%name, &
476 & 'lat', tout, &
477 & (/1,flt(ng)%Rindex/), &
478 & (/nfloats(ng),1/), &
479 & piofile = flt(ng)%pioFile, &
480 & piovar = flt(ng)%pioVar(idglat)%vd)
481 ELSE
482 CALL pio_netcdf_put_fvar (ng, inlm, flt(ng)%name, &
483 & 'y', tout, &
484 & (/1,flt(ng)%Rindex/), &
485 & (/nfloats(ng),1/), &
486 & piofile = flt(ng)%pioFile, &
487 & piovar = flt(ng)%pioVar(idglat)%vd)
488 END IF
489 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
490
491# ifdef SOLVE3D
492!
493! Write out floats depths.
494!
495 DO l=1,nfloats(ng)
496 IF (drifter(ng)%bounded(l)) THEN
497 tout(l)=drifter(ng)%track(idpth,nf(ng),l)
498 ELSE
499 tout(l)=spval
500 END IF
501 END DO
502 CALL pio_netcdf_put_fvar (ng, inlm, flt(ng)%name, &
503 & 'depth', tout, &
504 & (/1,flt(ng)%Rindex/), &
505 & (/nfloats(ng),1/), &
506 & piofile = flt(ng)%pioFile, &
507 & piovar = flt(ng)%pioVar(iddpth)%vd)
508 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
509!
510! Write out density anomaly.
511!
512 DO l=1,nfloats(ng)
513 IF (drifter(ng)%bounded(l)) THEN
514 tout(l)=drifter(ng)%track(ifden,nf(ng),l)
515 ELSE
516 tout(l)=spval
517 END IF
518 END DO
519 CALL pio_netcdf_put_fvar (ng, inlm, flt(ng)%name, &
520 & trim(vname(1,iddano)), tout, &
521 & (/1,flt(ng)%Rindex/), &
522 & (/nfloats(ng),1/), &
523 & piofile = flt(ng)%pioFile, &
524 & piovar = flt(ng)%pioVar(iddano)%vd)
525 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
526!
527! Write out tracer type variables.
528!
529 DO itrc=1,nt(ng)
530 DO l=1,nfloats(ng)
531 IF (drifter(ng)%bounded(l)) THEN
532 tout(l)=drifter(ng)%track(iftvar(itrc),nf(ng),l)
533 ELSE
534 tout(l)=spval
535 END IF
536 END DO
537 CALL pio_netcdf_put_fvar (ng, inlm, flt(ng)%name, &
538 & trim(vname(1,idtvar(itrc))), tout, &
539 & (/1,flt(ng)%Rindex/), &
540 & (/nfloats(ng),1/), &
541 & piofile = flt(ng)%pioFile, &
542 & piovar = flt(ng)%pioTrc(itrc)%vd)
543 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
544 END DO
545# endif
546
547# ifdef FLOAT_OYSTER
548!
549! Write out biological float swimming time.
550!
551 DO l=1,nfloats(ng)
552 IF (drifter(ng)%bounded(l)) THEN
553 tout(l)=drifter(ng)%track(iswim,nf(ng),l)
554 ELSE
555 tout(l)=spval
556 END IF
557 END DO
558 CALL pio_netcdf_put_fvar (ng, inlm, flt(ng)%name, &
559 & 'swim_time', tout, &
560 & (/1,flt(ng)%Rindex/), &
561 & (/nfloats(ng),1/), &
562 & piofile = flt(ng)%pioFile, &
563 & piovar = flt(ng)%pioVar(idswim)%vd)
564 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
565!
566! Write out biological float vertical velocity.
567!
568 DO l=1,nfloats(ng)
569 IF (drifter(ng)%bounded(l)) THEN
570 tout(l)=drifter(ng)%track(iwbio,nf(ng),l)
571 ELSE
572 tout(l)=spval
573 END IF
574 END DO
575 CALL pio_netcdf_put_fvar (ng, inlm, flt(ng)%name, &
576 & 'w_bio', tout, &
577 & (/1,flt(ng)%Rindex/), &
578 & (/nfloats(ng),1/), &
579 & piofile = flt(ng)%pioFile, &
580 & piovar = flt(ng)%pioVar(idwbio)%vd)
581 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
582!
583! Write out biological float size (length).
584!
585 DO l=1,nfloats(ng)
586 IF (drifter(ng)%bounded(l)) THEN
587 tout(l)=drifter(ng)%track(isizf,nf(ng),l)
588 ELSE
589 tout(l)=spval
590 END IF
591 END DO
592 CALL pio_netcdf_put_fvar (ng, inlm, flt(ng)%name, &
593 & 'bio_size', tout, &
594 & (/1,flt(ng)%Rindex/), &
595 & (/nfloats(ng),1/), &
596 & piofile = flt(ng)%pioFile, &
597 & piovar = flt(ng)%pioVar(idsize)%vd)
598 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
599!
600! Write out biological float sinking velocity.
601!
602 DO l=1,nfloats(ng)
603 IF (drifter(ng)%bounded(l)) THEN
604 tout(l)=drifter(ng)%track(iwsin,nf(ng),l)
605 ELSE
606 tout(l)=spval
607 END IF
608 END DO
609 CALL pio_netcdf_put_fvar (ng, inlm, flt(ng)%name, &
610 & 'bio_sink', tout, &
611 & (/1,flt(ng)%Rindex/), &
612 & (/nfloats(ng),1/), &
613 & piofile = flt(ng)%pioFile, &
614 & piovar = flt(ng)%pioVar(idwsin)%vd)
615 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
616# endif
617!
618!-----------------------------------------------------------------------
619! Synchronize floats NetCDF file to disk.
620!-----------------------------------------------------------------------
621!
622 CALL pio_netcdf_sync (ng, inlm, flt(ng)%name, flt(ng)%pioFile)
623!
624 RETURN
subroutine, public pio_netcdf_sync(ng, model, ncname, piofile)

References mod_floats::drifter, mod_scalars::exit_flag, mod_iounits::flt, strings_mod::founderror(), mod_ncparam::iddano, mod_ncparam::iddpth, mod_ncparam::idglat, mod_ncparam::idglon, mod_floats::idpth, mod_ncparam::idsize, mod_ncparam::idswim, mod_ncparam::idtime, mod_ncparam::idtvar, mod_ncparam::idwbio, mod_ncparam::idwsin, mod_ncparam::idxgrd, mod_ncparam::idygrd, mod_ncparam::idzgrd, mod_floats::ifden, mod_floats::iflat, mod_floats::iflon, mod_floats::iftvar, mod_param::inlm, mod_floats::isizf, mod_floats::iswim, mod_floats::iwbio, mod_floats::iwsin, mod_floats::ixgrd, mod_floats::iygrd, mod_floats::izgrd, mod_stepping::nf, mod_param::nfloats, mod_scalars::noerror, mod_param::nt, mod_pio_netcdf::pio_netcdf_sync(), mod_iounits::sourcefile, mod_scalars::spherical, mod_scalars::spval, mod_scalars::time, and mod_ncparam::vname.

Referenced by wrt_floats().

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