ROMS
Loading...
Searching...
No Matches
output.F
Go to the documentation of this file.
1#include "cppdefs.h"
2#ifdef NONLINEAR
3 SUBROUTINE output (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 subroutine manages nonlinear model output. It creates output !
13! NetCDF files and writes out data into NetCDF files. If requested, !
14! it can create several history and/or time-averaged files to avoid !
15! generating too large files during a single model run. !
16! !
17!=======================================================================
18!
19 USE mod_param
20 USE mod_parallel
21# ifdef FLOATS
22 USE mod_floats
23# endif
24# if defined FOUR_DVAR || defined VERIFICATION
25 USE mod_fourdvar
26# endif
27 USE mod_iounits
28 USE mod_ncparam
29 USE mod_scalars
30!
31 USE close_io_mod, ONLY : close_file
32# ifdef AVERAGES
33 USE def_avg_mod, ONLY : def_avg
34# endif
35# ifdef DIAGNOSTICS
36 USE def_diags_mod, ONLY : def_diags
37# endif
38# ifdef GRID_EXTRACT
39 USE def_extract_mod, ONLY : def_extract
40# endif
41# ifdef FLOATS
42 USE def_floats_mod, ONLY : def_floats
43# endif
44 USE def_his_mod, ONLY : def_his
45 USE def_quick_mod, ONLY : def_quick
46 USE def_rst_mod, ONLY : def_rst
47# ifdef STATIONS
49# endif
50# ifdef DISTRIBUTE
51 USE distribute_mod, ONLY : mp_bcasts
52# endif
53# ifdef OBSERVATIONS
54 USE obs_read_mod, ONLY : obs_read
55 USE obs_write_mod, ONLY : obs_write
56# endif
57 USE strings_mod, ONLY : founderror
58# ifdef AVERAGES
59 USE wrt_avg_mod, ONLY : wrt_avg
60# endif
61# ifdef DIAGNOSTICS
62 USE wrt_diags_mod, ONLY : wrt_diags
63# endif
64# ifdef GRID_EXTRACT
65 USE wrt_extract_mod, ONLY : wrt_extract
66# endif
67# ifdef FLOATS
68 USE wrt_floats_mod, ONLY : wrt_floats
69# endif
70 USE wrt_his_mod, ONLY : wrt_his
71 USE wrt_quick_mod, ONLY : wrt_quick
72 USE wrt_rst_mod, ONLY : wrt_rst
73# ifdef STATIONS
75# endif
76# if defined AVERAGES && defined AVERAGES_DETIDE && \
77 (defined ssh_tides || defined uv_tides)
78 USE wrt_tides_mod, ONLY : wrt_tides
79# endif
80!
81 implicit none
82!
83! Imported variable declarations.
84!
85 integer, intent(in) :: ng
86!
87! Local variable declarations.
88!
89 logical :: Ldefine, Lupdate, NewFile
90!
91 integer :: Fcount, ifile, status, tile
92!
93 character (len=*), parameter :: MyFile = &
94 & __FILE__
95!
96 sourcefile=myfile
97
98# ifdef PROFILE
99!
100!-----------------------------------------------------------------------
101! Turn on output data time wall clock.
102!-----------------------------------------------------------------------
103!
104 CALL wclock_on (ng, inlm, 8, __line__, myfile)
105# endif
106!
107!-----------------------------------------------------------------------
108! If appropriate, process nonlinear history NetCDF file.
109!-----------------------------------------------------------------------
110!
111! Set tile for local array manipulations in output routines.
112!
113# ifdef DISTRIBUTE
114 tile=myrank
115# else
116 tile=-1
117# endif
118!
119! Turn off checking for analytical header files.
120!
121 IF (lanafile) THEN
122 lanafile=.false.
123 END IF
124!
125! If appropriate, set switch for updating biology header file global
126! attribute in output NetCDF files.
127!
128# ifdef BIOLOGY
129 lupdate=.true.
130# else
131 lupdate=.false.
132# endif
133!
134! Create output history NetCDF file or prepare existing file to
135! append new data to it. Also, notice that it is possible to
136! create several files during a single model run.
137!
138 IF (ldefhis(ng)) THEN
139 IF (ndefhis(ng).gt.0) THEN
140 IF (idefhis(ng).lt.0) THEN
141 idefhis(ng)=((ntstart(ng)-1)/ndefhis(ng))*ndefhis(ng)
142 IF (idefhis(ng).lt.iic(ng)-1) THEN
143 idefhis(ng)=idefhis(ng)+ndefhis(ng)
144 END IF
145 END IF
146 IF ((nrrec(ng).ne.0).and.(iic(ng).eq.ntstart(ng))) THEN
147 IF ((iic(ng)-1).eq.idefhis(ng)) THEN
148 his(ng)%load=0 ! restart, reset counter
149 ldefine=.false. ! finished file, delay
150 ELSE ! creation of next file
151 ldefine=.true.
152 newfile=.false. ! unfinished file, inquire
153 END IF ! content for appending
154 idefhis(ng)=idefhis(ng)+nhis(ng) ! restart offset
155 ELSE IF ((iic(ng)-1).eq.idefhis(ng)) THEN
156 idefhis(ng)=idefhis(ng)+ndefhis(ng)
157 IF (nhis(ng).ne.ndefhis(ng).and.iic(ng).eq.ntstart(ng)) THEN
158 idefhis(ng)=idefhis(ng)+nhis(ng) ! multiple record offset
159 END IF
160 ldefine=.true.
161 newfile=.true.
162 ELSE
163 ldefine=.false.
164 END IF
165 IF (ldefine) THEN ! create new file or
166 IF (iic(ng).eq.ntstart(ng)) THEN ! inquire existing file
167 his(ng)%load=0 ! reset filename counter
168 END IF
169 ifile=(iic(ng)-1)/ndefhis(ng)+1 ! next filename suffix
170 his(ng)%load=his(ng)%load+1
171 IF (his(ng)%load.gt.his(ng)%Nfiles) THEN
172 IF (master) THEN
173 WRITE (stdout,10) 'HIS(ng)%load = ', his(ng)%load, &
174 & his(ng)%Nfiles, trim(his(ng)%base), &
175 & ifile
176 END IF
177 exit_flag=4
179 & __line__, myfile)) RETURN
180 END IF
181 fcount=his(ng)%load
182 his(ng)%Nrec(fcount)=0
183 IF (master) THEN
184 WRITE (his(ng)%name,20) trim(his(ng)%base), ifile
185 END IF
186# ifdef DISTRIBUTE
187 CALL mp_bcasts (ng, inlm, his(ng)%name)
188# endif
189 his(ng)%files(fcount)=trim(his(ng)%name)
190 CALL close_file (ng, inlm, his(ng), his(ng)%name, lupdate)
191 CALL def_his (ng, newfile)
192 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
193 END IF
194 IF ((iic(ng).eq.ntstart(ng)).and.(nrrec(ng).ne.0)) THEN
195 lwrthis(ng)=.false. ! avoid writing initial
196 ELSE ! fields during restart
197 lwrthis(ng)=.true.
198 END IF
199 ELSE
200 IF (iic(ng).eq.ntstart(ng)) THEN
201 CALL def_his (ng, ldefout(ng))
202 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
203 lwrthis(ng)=.true.
204 ldefhis(ng)=.false.
205 END IF
206 END IF
207 END IF
208!
209! Write out data into history NetCDF file. Avoid writing initial
210! conditions in perturbation mode computations.
211!
212 IF (lwrthis(ng)) THEN
213 IF (lwrtper(ng)) THEN
214 IF ((iic(ng).gt.ntstart(ng)).and. &
215 & (mod(iic(ng)-1,nhis(ng)).eq.0)) THEN
216 IF (nrrec(ng).eq.0.or.iic(ng).ne.ntstart(ng)) THEN
217 CALL wrt_his (ng, tile)
218 END IF
219 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
220 END IF
221 ELSE
222 IF (mod(iic(ng)-1,nhis(ng)).eq.0) THEN
223 CALL wrt_his (ng, tile)
224 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
225 END IF
226 END IF
227 END IF
228
229# ifdef GRID_EXTRACT
230!
231!-----------------------------------------------------------------------
232! If appropriate, process nonlinear history extraction NetCDF file.
233!-----------------------------------------------------------------------
234!
235! Create output history extraction NetCDF file or prepare existing
236! file to append new data to it. Also, notice that it is possible
237! to create several files during a single model run.
238!
239 IF (ldefxtr(ng)) THEN
240 IF (ndefxtr(ng).gt.0) THEN
241 IF (idefxtr(ng).lt.0) THEN
242 idefxtr(ng)=((ntstart(ng)-1)/ndefxtr(ng))*ndefxtr(ng)
243 IF (idefxtr(ng).lt.iic(ng)-1) THEN
244 idefxtr(ng)=idefxtr(ng)+ndefxtr(ng)
245 END IF
246 END IF
247 IF ((nrrec(ng).ne.0).and.(iic(ng).eq.ntstart(ng))) THEN
248 IF ((iic(ng)-1).eq.idefxtr(ng)) THEN
249 xtr(ng)%load=0 ! restart, reset counter
250 ldefine=.false. ! finished file, delay
251 ELSE ! creation of next file
252 ldefine=.true.
253 newfile=.false. ! unfinished file, inquire
254 END IF ! content for appending
255 idefxtr(ng)=idefxtr(ng)+nxtr(ng) ! restart offset
256 ELSE IF ((iic(ng)-1).eq.idefxtr(ng)) THEN
257 idefxtr(ng)=idefxtr(ng)+ndefxtr(ng)
258 IF (nxtr(ng).ne.ndefxtr(ng).and.iic(ng).eq.ntstart(ng)) THEN
259 idefxtr(ng)=idefxtr(ng)+nxtr(ng) ! multiple record offset
260 END IF
261 ldefine=.true.
262 newfile=.true.
263 ELSE
264 ldefine=.false.
265 END IF
266 IF (ldefine) THEN ! create new file or
267 IF (iic(ng).eq.ntstart(ng)) THEN ! inquire existing file
268 xtr(ng)%load=0 ! reset filename counter
269 END IF
270 ifile=(iic(ng)-1)/ndefxtr(ng)+1 ! next filename suffix
271 xtr(ng)%load=xtr(ng)%load+1
272 IF (xtr(ng)%load.gt.xtr(ng)%Nfiles) THEN
273 IF (master) THEN
274 WRITE (stdout,10) 'XTR(ng)%load = ', xtr(ng)%load, &
275 & xtr(ng)%Nfiles, trim(xtr(ng)%base), &
276 & ifile
277 END IF
278 exit_flag=4
280 & __line__, myfile)) RETURN
281 END IF
282 fcount=xtr(ng)%load
283 xtr(ng)%Nrec(fcount)=0
284 IF (master) THEN
285 WRITE (xtr(ng)%name,20) trim(xtr(ng)%base), ifile
286 END IF
287# ifdef DISTRIBUTE
288 CALL mp_bcasts (ng, inlm, xtr(ng)%name)
289# endif
290 xtr(ng)%files(fcount)=trim(xtr(ng)%name)
291 CALL close_file (ng, inlm, xtr(ng), xtr(ng)%name, lupdate)
292 CALL def_extract (ng, newfile)
293 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
294 END IF
295 IF ((iic(ng).eq.ntstart(ng)).and.(nrrec(ng).ne.0)) THEN
296 lwrtxtr(ng)=.false. ! avoid writing initial
297 ELSE ! fields during restart
298 lwrtxtr(ng)=.true.
299 END IF
300 ELSE
301 IF (iic(ng).eq.ntstart(ng)) THEN
302 CALL def_extract (ng, ldefout(ng))
303 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
304 lwrtxtr(ng)=.true.
305 ldefxtr(ng)=.false.
306 END IF
307 END IF
308 END IF
309!
310! Write out data into history NetCDF file. Avoid writing initial
311! conditions in perturbation mode computations.
312!
313 IF (lwrtxtr(ng)) THEN
314 IF (lwrtper(ng)) THEN
315 IF ((iic(ng).gt.ntstart(ng)).and. &
316 & (mod(iic(ng)-1,nxtr(ng)).eq.0)) THEN
317 IF (nrrec(ng).eq.0.or.iic(ng).ne.ntstart(ng)) THEN
318 CALL wrt_extract (ng, tile)
319 END IF
320 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
321 END IF
322 ELSE
323 IF (mod(iic(ng)-1,nxtr(ng)).eq.0) THEN
324 CALL wrt_extract (ng, tile)
325 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
326 END IF
327 END IF
328 END IF
329# endif
330!
331!-----------------------------------------------------------------------
332! If appropriate, process nonlinear quicksave NetCDF file.
333!-----------------------------------------------------------------------
334!
335! Create output quicksave NetCDF file or prepare existing file to
336! append new data to it. Also, notice that it is possible to
337! create several files during a single model run.
338!
339 IF (ldefqck(ng)) THEN
340 IF (ndefqck(ng).gt.0) THEN
341 IF (idefqck(ng).lt.0) THEN
342 idefqck(ng)=((ntstart(ng)-1)/ndefqck(ng))*ndefqck(ng)
343 IF (idefqck(ng).lt.iic(ng)-1) THEN
344 idefqck(ng)=idefqck(ng)+ndefqck(ng)
345 END IF
346 END IF
347 IF ((nrrec(ng).ne.0).and.(iic(ng).eq.ntstart(ng))) THEN
348 IF ((iic(ng)-1).eq.idefqck(ng)) THEN
349 qck(ng)%load=0 ! restart, reset counter
350 ldefine=.false. ! finished file, delay
351 ELSE ! creation of next file
352 ldefine=.true.
353 newfile=.false. ! unfinished file, inquire
354 END IF ! content for appending
355 idefqck(ng)=idefqck(ng)+nqck(ng) ! restart offset
356 ELSE IF ((iic(ng)-1).eq.idefqck(ng)) THEN
357 idefqck(ng)=idefqck(ng)+ndefqck(ng)
358 IF (nqck(ng).ne.ndefqck(ng).and.iic(ng).eq.ntstart(ng)) THEN
359 idefqck(ng)=idefqck(ng)+nqck(ng) ! multiple record offset
360 END IF
361 ldefine=.true.
362 newfile=.true.
363 ELSE
364 ldefine=.false.
365 END IF
366 IF (ldefine) THEN ! create new file or
367 IF (iic(ng).eq.ntstart(ng)) THEN ! inquire existing file
368 qck(ng)%load=0 ! reset filename counter
369 END IF
370 ifile=(iic(ng)-1)/ndefqck(ng)+1 ! next filename suffix
371 qck(ng)%load=qck(ng)%load+1
372 IF (qck(ng)%load.gt.qck(ng)%Nfiles) THEN
373 IF (master) THEN
374 WRITE (stdout,10) 'QCK(ng)%load = ', qck(ng)%load, &
375 & qck(ng)%Nfiles, trim(qck(ng)%base), &
376 & ifile
377 END IF
378 exit_flag=4
380 & __line__, myfile)) RETURN
381 END IF
382 fcount=qck(ng)%load
383 qck(ng)%Nrec(fcount)=0
384 IF (master) THEN
385 WRITE (qck(ng)%name,20) trim(qck(ng)%base), ifile
386 END IF
387# ifdef DISTRIBUTE
388 CALL mp_bcasts (ng, inlm, qck(ng)%name)
389# endif
390 qck(ng)%files(fcount)=trim(qck(ng)%name)
391 CALL close_file (ng, inlm, qck(ng), qck(ng)%name, lupdate)
392 CALL def_quick (ng, newfile)
393 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
394 END IF
395 IF ((iic(ng).eq.ntstart(ng)).and.(nrrec(ng).ne.0)) THEN
396 lwrtqck(ng)=.false. ! avoid writing initial
397 ELSE ! fields during restart
398 lwrtqck(ng)=.true.
399 END IF
400 ELSE
401 IF (iic(ng).eq.ntstart(ng)) THEN
402 CALL def_quick (ng, ldefout(ng))
403 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
404 lwrtqck(ng)=.true.
405 ldefqck(ng)=.false.
406 END IF
407 END IF
408 END IF
409!
410! Write out data into quicksave NetCDF file. Avoid writing initial
411! conditions in perturbation mode computations.
412!
413 IF (lwrtqck(ng)) THEN
414 IF (lwrtper(ng)) THEN
415 IF ((iic(ng).gt.ntstart(ng)).and. &
416 & (mod(iic(ng)-1,nqck(ng)).eq.0)) THEN
417 IF (nrrec(ng).eq.0.or.iic(ng).ne.ntstart(ng)) THEN
418 CALL wrt_quick (ng, tile)
419 END IF
420 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
421 END IF
422 ELSE
423 IF (mod(iic(ng)-1,nqck(ng)).eq.0) THEN
424 CALL wrt_quick (ng, tile)
425 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
426 END IF
427 END IF
428 END IF
429
430# ifdef AVERAGES
431!
432!-----------------------------------------------------------------------
433! If appropriate, process time-averaged NetCDF file.
434!-----------------------------------------------------------------------
435!
436! Create output time-averaged NetCDF file or prepare existing file
437! to append new data to it. Also, notice that it is possible to
438! create several files during a single model run.
439!
440 IF (ldefavg(ng)) THEN
441 IF (ndefavg(ng).gt.0) THEN
442 IF (idefavg(ng).lt.0) THEN
443 idefavg(ng)=((ntstart(ng)-1)/ndefavg(ng))*ndefavg(ng)
444 IF ((ndefavg(ng).eq.navg(ng)).and.(idefavg(ng).le.0)) THEN
445 idefavg(ng)=ndefavg(ng) ! one file per record
446 ELSE IF (idefavg(ng).lt.iic(ng)-1) THEN
447 idefavg(ng)=idefavg(ng)+ndefavg(ng)
448 END IF
449 END IF
450 IF ((nrrec(ng).ne.0).and.(iic(ng).eq.ntstart(ng))) THEN
451 IF ((iic(ng)-1).eq.idefavg(ng)) THEN
452 avg(ng)%load=0 ! restart, reset counter
453 ldefine=.false. ! finished file, delay
454 ELSE ! creation of next file
455 newfile=.false.
456 ldefine=.true. ! unfinished file, inquire
457 END IF ! content for appending
458 idefavg(ng)=idefavg(ng)+navg(ng) ! restart offset
459 ELSE IF ((iic(ng)-1).eq.idefavg(ng)) THEN
460 idefavg(ng)=idefavg(ng)+ndefavg(ng)
461 IF (navg(ng).ne.ndefavg(ng).and.iic(ng).eq.ntstart(ng)) THEN
462 idefavg(ng)=idefavg(ng)+navg(ng)
463 END IF
464 ldefine=.true.
465 newfile=.true.
466 ELSE
467 ldefine=.false.
468 END IF
469 IF (ldefine) THEN
470 IF (iic(ng).eq.ntstart(ng)) THEN
471 avg(ng)%load=0 ! reset filename counter
472 END IF
473 IF (ndefavg(ng).eq.navg(ng)) THEN ! next filename suffix
474 ifile=(iic(ng)-1)/ndefavg(ng)
475 ELSE
476 ifile=(iic(ng)-1)/ndefavg(ng)+1
477 END IF
478 avg(ng)%load=avg(ng)%load+1
479 IF (avg(ng)%load.gt.avg(ng)%Nfiles) THEN
480 IF (master) THEN
481 WRITE (stdout,10) 'AVG(ng)%load = ', avg(ng)%load, &
482 & avg(ng)%Nfiles, trim(avg(ng)%base), &
483 & ifile
484 END IF
485 exit_flag=4
487 & __line__, myfile)) RETURN
488 END IF
489 fcount=avg(ng)%load
490 avg(ng)%Nrec(fcount)=0
491 IF (master) THEN
492 WRITE (avg(ng)%name,20) trim(avg(ng)%base), ifile
493 END IF
494# ifdef DISTRIBUTE
495 CALL mp_bcasts (ng, inlm, avg(ng)%name)
496# endif
497 avg(ng)%files(fcount)=trim(avg(ng)%name)
498 CALL close_file (ng, inlm, avg(ng), avg(ng)%name, lupdate)
499 CALL def_avg (ng, newfile)
500 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
501 lwrtavg(ng)=.true.
502 END IF
503 ELSE
504 IF (iic(ng).eq.ntstart(ng)) THEN
505 CALL def_avg (ng, ldefout(ng))
506 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
507 lwrtavg(ng)=.true.
508 ldefavg(ng)=.false.
509 END IF
510 END IF
511 END IF
512!
513! Write out data into time-averaged NetCDF file.
514!
515 IF (lwrtavg(ng)) THEN
516 IF (((iic(ng).gt.ntstart(ng)).and. &
517 & (mod(iic(ng)-1,navg(ng)).eq.0)).or. &
518 & ((iic(ng).ge.ntsavg(ng)).and.(navg(ng).eq.1))) THEN
519 CALL wrt_avg (ng, tile)
520 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
521# if defined AVERAGES_DETIDE && (defined SSH_TIDES || defined UV_TIDES)
522 CALL wrt_tides (ng, tile)
523 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
524# endif
525 END IF
526 END IF
527# endif
528
529# ifdef DIAGNOSTICS
530!
531!-----------------------------------------------------------------------
532! If appropriate, process time-averaged diagnostics NetCDF file.
533!-----------------------------------------------------------------------
534!
535! Create output time-averaged diagnostics NetCDF file or prepare
536! existing file to append new data to it. Also, notice that it is
537! possible to create several files during a single model run.
538!
539 IF (ldefdia(ng)) THEN
540 IF (ndefdia(ng).gt.0) THEN
541 IF (idefdia(ng).lt.0) THEN
542 idefdia(ng)=((ntstart(ng)-1)/ndefdia(ng))*ndefdia(ng)
543 IF ((ndefdia(ng).eq.ndia(ng)).and.(idefdia(ng).le.0)) THEN
544 idefdia(ng)=ndefdia(ng) ! one file per record
545 ELSE IF (idefdia(ng).lt.iic(ng)-1) THEN
546 idefdia(ng)=idefdia(ng)+ndefdia(ng)
547 END IF
548 END IF
549 IF ((nrrec(ng).ne.0).and.(iic(ng).eq.ntstart(ng))) THEN
550 IF ((iic(ng)-1).eq.idefdia(ng)) THEN
551 dia(ng)%load=0 ! restart, reset counter
552 ldefine=.false. ! finished file, delay
553 ELSE ! creation of next file
554 newfile=.false.
555 ldefine=.true. ! unfinished file, inquire
556 END IF ! content for appending
557 idefdia(ng)=idefdia(ng)+ndia(ng) ! restart offset
558 ELSE IF ((iic(ng)-1).eq.idefdia(ng)) THEN
559 idefdia(ng)=idefdia(ng)+ndefdia(ng)
560 IF (ndia(ng).ne.ndefdia(ng).and.iic(ng).eq.ntstart(ng)) THEN
561 idefdia(ng)=idefdia(ng)+ndia(ng)
562 END IF
563 ldefine=.true.
564 newfile=.true.
565 ELSE
566 ldefine=.false.
567 END IF
568 IF (ldefine) THEN
569 IF (iic(ng).eq.ntstart(ng)) THEN
570 dia(ng)%load=0 ! reset filename counter
571 END IF
572 IF (ndefdia(ng).eq.ndia(ng)) THEN ! next filename suffix
573 ifile=(iic(ng)-1)/ndefdia(ng)
574 ELSE
575 ifile=(iic(ng)-1)/ndefdia(ng)+1
576 END IF
577 dia(ng)%load=dia(ng)%load+1
578 IF (dia(ng)%load.gt.dia(ng)%Nfiles) THEN
579 IF (master) THEN
580 WRITE (stdout,10) 'DIA(ng)%load = ', dia(ng)%load, &
581 & dia(ng)%Nfiles, trim(dia(ng)%base), &
582 & ifile
583 END IF
584 exit_flag=4
586 & __line__, myfile)) RETURN
587 END IF
588 fcount=dia(ng)%load
589 dia(ng)%Nrec(fcount)=0
590 IF (master) THEN
591 WRITE (dia(ng)%name,20) trim(dia(ng)%base), ifile
592 END IF
593# ifdef DISTRIBUTE
594 CALL mp_bcasts (ng, inlm, dia(ng)%name)
595# endif
596 dia(ng)%files(fcount)=trim(dia(ng)%name)
597 CALL close_file (ng, inlm, dia(ng), dia(ng)%name, lupdate)
598 CALL def_diags (ng, newfile)
599 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
600 lwrtdia(ng)=.true.
601 END IF
602 ELSE
603 IF (iic(ng).eq.ntstart(ng)) THEN
604 CALL def_diags (ng, ldefout(ng))
605 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
606 lwrtdia(ng)=.true.
607 ldefdia(ng)=.false.
608 END IF
609 END IF
610 END IF
611!
612! Write out data into time-averaged diagnostics NetCDF file.
613!
614 IF (lwrtdia(ng)) THEN
615 IF (((iic(ng).gt.ntstart(ng)).and. &
616 & (mod(iic(ng)-1,ndia(ng)).eq.0)).or. &
617 & ((iic(ng).ge.ntsdia(ng)).and.(ndia(ng).eq.1))) THEN
618 CALL wrt_diags (ng, tile)
619 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
620 END IF
621 END IF
622# endif
623
624# ifdef STATIONS
625!
626!-----------------------------------------------------------------------
627! If appropriate, process stations NetCDF file.
628!-----------------------------------------------------------------------
629!
630 IF (lstations(ng).and. &
631 & (nstation(ng).gt.0).and.(nsta(ng).gt.0)) THEN
632!
633! Create output station NetCDF file or prepare existing file to
634! append new data to it.
635!
636 IF (ldefsta(ng).and.(iic(ng).eq.ntstart(ng))) THEN
637 CALL def_station (ng, ldefout(ng))
638 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
639 ldefsta(ng)=.false.
640 END IF
641!
642! Write out data into stations NetCDF file.
643!
644 IF (mod(iic(ng)-1,nsta(ng)).eq.0) THEN
645 CALL wrt_station (ng, tile)
646 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
647 END IF
648 END IF
649# endif
650
651# ifdef FLOATS
652!
653!-----------------------------------------------------------------------
654! If appropriate, process floats NetCDF file.
655!-----------------------------------------------------------------------
656!
657 IF (lfloats(ng).and. &
658 & (nfloats(ng).gt.0).and.(nflt(ng).gt.0)) THEN
659!
660! Create output floats NetCDF file or prepare existing file to
661! append new data to it.
662!
663 IF (ldefflt(ng)) THEN
664 IF (frrec(ng).eq.0) THEN
665 newfile=.true.
666 ELSE
667 newfile=.false.
668 END IF
669 CALL def_floats (ng, newfile)
670 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
671 ldefflt(ng)=.false.
672 END IF
673!
674! Write out data into floats NetCDF file.
675!
676 IF ((mod(iic(ng)-1,nflt(ng)).eq.0).and. &
677 & ((frrec(ng).eq.0).or.(iic(ng).ne.ntstart(ng)))) THEN
678 CALL wrt_floats (ng)
679 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
680 END IF
681 END IF
682# endif
683!
684!-----------------------------------------------------------------------
685! If appropriate, process restart NetCDF file.
686!-----------------------------------------------------------------------
687!
688! Create output restart NetCDF file or prepare existing file to
689! append new data to it.
690!
691 IF (ldefrst(ng)) THEN
692 CALL def_rst (ng)
693 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
694 lwrtrst(ng)=.true.
695 ldefrst(ng)=.false.
696 END IF
697!
698! Write out data into restart NetCDF file.
699!
700 IF (lwrtrst(ng)) THEN
701 IF ((iic(ng).gt.ntstart(ng)).and. &
702 & (mod(iic(ng)-1,nrst(ng)).eq.0)) THEN
703 CALL wrt_rst (ng, tile)
704 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
705 END IF
706 END IF
707
708# if (defined FOUR_DVAR || \
709 defined verification) && \
710 !defined I4DVAR_ANA_SENSITIVITY
711# ifdef OBSERVATIONS
712!
713!-----------------------------------------------------------------------
714! If appropriate, process and write model state at observation
715! locations. Compute misfit (model-observations) cost function.
716!-----------------------------------------------------------------------
717!
718 IF (((time(ng)-0.5_r8*dt(ng)).le.obstime(ng)).and. &
719 & (obstime(ng).lt.(time(ng)+0.5_r8*dt(ng)))) THEN
720 processobs=.true.
721 CALL obs_read (ng, inlm, .false.)
722# ifdef SP4DVAR
723 IF (iic(ng).ne.ntend(ng)+1.or.lsadd(ng)) THEN
724 CALL obs_write (ng, tile, inlm)
725 CALL obs_cost (ng, inlm)
726 END IF
727# else
728 CALL obs_write (ng, tile, inlm)
729# if !(defined R4DVAR || defined VERIFICATION)
730 CALL obs_cost (ng, inlm)
731# endif
732# endif
733 ELSE
734 processobs=.false.
735 END IF
736# endif
737# endif
738# ifdef PROFILE
739!
740!-----------------------------------------------------------------------
741! Turn off output data time wall clock.
742!-----------------------------------------------------------------------
743!
744 CALL wclock_off (ng, inlm, 8, __line__, myfile)
745# endif
746!
747 10 FORMAT (/,' OUTPUT - multi-file counter ',a,i0, &
748 & ', is greater than Nfiles = ',i0,1x,'dimension', &
749 & /,10x,'in structure when creating next file: ', &
750 & a,'_',i4.4,'.nc', &
751 & /,10x,'Incorrect OutFiles logic in ''read_phypar''.')
752 20 FORMAT (a,'_',i4.4,'.nc')
753!
754 RETURN
755 END SUBROUTINE output
756#else
757 SUBROUTINE output
758 RETURN
759 END SUBROUTINE output
760#endif
subroutine, public close_file(ng, model, s, ncname, lupdate)
Definition close_io.F:43
subroutine, public def_avg(ng, ldef)
Definition def_avg.F:79
subroutine, public def_diags(ng, ldef)
Definition def_diags.F:51
subroutine, public def_floats(ng, ldef)
Definition def_floats.F:53
subroutine, public def_his(ng, ldef)
Definition def_his.F:74
subroutine, public def_quick(ng, ldef)
Definition def_quick.F:74
subroutine, public def_rst(ng)
Definition def_rst.F:61
subroutine, public def_station(ng, ldef)
Definition def_station.F:75
integer, dimension(:), allocatable frrec
Definition mod_floats.F:138
logical, dimension(:), allocatable lsadd
logical, dimension(:), allocatable processobs
type(t_io), dimension(:), allocatable his
type(t_io), dimension(:), allocatable xtr
type(t_io), dimension(:), allocatable qck
type(t_io), dimension(:), allocatable avg
integer stdout
character(len=256) sourcefile
type(t_io), dimension(:), allocatable dia
integer, dimension(:), allocatable idefqck
integer, dimension(:), allocatable idefavg
integer, dimension(:), allocatable idefdia
logical lanafile
integer, dimension(:), allocatable idefxtr
integer, dimension(:), allocatable idefhis
logical master
integer, dimension(:), allocatable nfloats
Definition mod_param.F:543
integer, parameter inlm
Definition mod_param.F:662
integer, dimension(:), allocatable nstation
Definition mod_param.F:557
logical, dimension(:), allocatable lwrtqck
logical, dimension(:), allocatable lfloats
logical, dimension(:), allocatable lwrtdia
integer, dimension(:), allocatable nrrec
real(dp), dimension(:), allocatable obstime
integer, dimension(:), allocatable nxtr
integer, dimension(:), allocatable iic
integer, dimension(:), allocatable ndefhis
real(dp), dimension(:), allocatable dt
logical, dimension(:), allocatable lwrtxtr
logical, dimension(:), allocatable ldefflt
logical, dimension(:), allocatable ldefdia
integer, dimension(:), allocatable nrst
integer, dimension(:), allocatable nqck
logical, dimension(:), allocatable ldefavg
logical, dimension(:), allocatable ldefqck
integer, dimension(:), allocatable nflt
integer, dimension(:), allocatable navg
logical, dimension(:), allocatable lwrtavg
logical, dimension(:), allocatable ldefhis
integer, dimension(:), allocatable ntend
integer exit_flag
integer, dimension(:), allocatable ndefqck
integer, dimension(:), allocatable nsta
integer, dimension(:), allocatable nhis
integer, dimension(:), allocatable ndefavg
logical, dimension(:), allocatable lwrtper
logical, dimension(:), allocatable ldefxtr
integer, dimension(:), allocatable ndefxtr
logical, dimension(:), allocatable lwrthis
logical, dimension(:), allocatable ldefrst
logical, dimension(:), allocatable ldefout
logical, dimension(:), allocatable lwrtrst
real(dp), dimension(:), allocatable time
logical, dimension(:), allocatable ldefsta
logical, dimension(:), allocatable lstations
integer, dimension(:), allocatable ntstart
integer, dimension(:), allocatable ndia
integer, dimension(:), allocatable ntsdia
integer, dimension(:), allocatable ntsavg
integer noerror
integer, dimension(:), allocatable ndefdia
subroutine, public obs_read(ng, model, backward)
Definition obs_read.F:42
subroutine, public obs_write(ng, tile, model)
Definition obs_write.F:56
logical function, public founderror(flag, noerr, line, routine)
Definition strings.F:52
subroutine, public wrt_avg(ng, tile)
Definition wrt_avg.F:83
subroutine, public wrt_diags(ng, tile)
Definition wrt_diags.F:51
subroutine, public wrt_floats(ng)
Definition wrt_floats.F:40
subroutine, public wrt_his(ng, tile)
Definition wrt_his.F:98
subroutine, public wrt_quick(ng, tile)
Definition wrt_quick.F:88
subroutine, public wrt_rst(ng, tile)
Definition wrt_rst.F:63
subroutine, public wrt_station(ng, tile)
Definition wrt_station.F:81
subroutine, public wrt_tides(ng, tile)
Definition wrt_tides.F:46
subroutine obs_cost(ng, model)
Definition obs_cost.F:4
subroutine output(ng)
Definition output.F:4
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