ROMS
Loading...
Searching...
No Matches
get_idata.F File Reference
#include "cppdefs.h"
Include dependency graph for get_idata.F:

Go to the source code of this file.

Functions/Subroutines

subroutine get_idata (ng)
 

Function/Subroutine Documentation

◆ get_idata()

subroutine get_idata ( integer, intent(in) ng)

Definition at line 2 of file get_idata.F.

3!
4!git $Id$
5!================================================== Hernan G. Arango ===
6! Copyright (c) 2002-2025 The ROMS Group !
7! Licensed under a MIT/X style license !
8! See License_ROMS.md !
9!=======================================================================
10! !
11! This routine reads input data that needs to be obtained only once. !
12! !
13! Currently, this routine is only executed in serial mode by the !
14! main thread. !
15! !
16!=======================================================================
17!
18 USE mod_param
19 USE mod_parallel
20#ifdef RED_TIDE
21 USE mod_biology
22#endif
23 USE mod_grid
24 USE mod_iounits
25 USE mod_ncparam
26#if defined AVERAGES && defined AVERAGES_DETIDE && \
27 (defined ssh_tides || defined uv_tides)
28 USE mod_netcdf
29# if defined PIO_LIB && defined DISTRIBUTE
31# endif
32#endif
33#ifdef RED_TIDE
34 USE mod_ocean
35#endif
36 USE mod_scalars
37 USE mod_sources
38 USE mod_stepping
39#if defined SSH_TIDES || defined UV_TIDES
40 USE mod_tides
41#endif
42!
43#if defined AVERAGES && defined AVERAGES_DETIDE && \
44 (defined ssh_tides || defined uv_tides)
45 USE def_tides_mod, ONLY : def_tides
46#endif
47 USE nf_fread3d_mod, ONLY : nf_fread3d
48#ifdef SOLVE3D
49 USE nf_fread4d_mod, ONLY : nf_fread4d
50#endif
51 USE strings_mod, ONLY : founderror
52!
53 implicit none
54!
55! Imported variable declarations.
56!
57 integer, intent(in) :: ng
58!
59! Local variable declarations.
60!
61 logical, save :: recordless = .true.
62
63 logical, dimension(3) :: update = &
64 & (/ .FALSE., .FALSE., .FALSE. /)
65!
66 integer :: LBi, UBi, LBj, UBj
67 integer :: itrc, is
68
69#if defined AVERAGES && defined AVERAGES_DETIDE && \
70 (defined ssh_tides || defined uv_tides)
71 integer :: gtype, status, varid, Vsize(4)
72!
73 real(r8), parameter :: Fscl = 1.0_r8
74
75 real(r8) :: Fmin, Fmax, Htime
76#endif
77 real(r8) :: time_save = 0.0_r8
78!
79 character (len=*), parameter :: MyFile = &
80 & __FILE__
81
82#if defined AVERAGES && defined AVERAGES_DETIDE && \
83 (defined ssh_tides || defined uv_tides)
84# if defined PIO_LIB && defined DISTRIBUTE
85!
86 TYPE (IO_Desc_t), pointer :: ioDesc
87 TYPE (My_VarDesc), pointer :: pioVar
88# endif
89#endif
90!
91 sourcefile=myfile
92!
93! Lower and upper bounds for tiled arrays.
94!
95 lbi=lbound(grid(ng)%h,dim=1)
96 ubi=ubound(grid(ng)%h,dim=1)
97 lbj=lbound(grid(ng)%h,dim=2)
98 ubj=ubound(grid(ng)%h,dim=2)
99
100#if defined AVERAGES && defined AVERAGES_DETIDE && \
101 (defined ssh_tides || defined uv_tides)
102!
103! Set Vsize to zero to deactivate interpolation of input data to model
104! grid in "nf_fread2d" and "nf_fread3d".
105!
106 DO is=1,4
107 vsize(is)=0
108 END DO
109#endif
110
111#ifdef PROFILE
112!
113!-----------------------------------------------------------------------
114! Turn on input data time wall clock.
115!-----------------------------------------------------------------------
116!
117 CALL wclock_on (ng, inlm, 3, __line__, myfile)
118#endif
119
120#if defined SSH_TIDES || defined UV_TIDES
121!
122!-----------------------------------------------------------------------
123! Tide period, amplitude, phase, and currents.
124!-----------------------------------------------------------------------
125!
126! Tidal Period.
127!
128 IF (lprocesstides(ng)) THEN
129 IF (iic(ng).eq.0) THEN
130 CALL get_ngfld (ng, inlm, idtper, tide(ng)%ncid, &
131# if defined PIO_LIB && defined DISTRIBUTE
132 & tide(ng)%pioFile, &
133# endif
134 & 1, tide(ng), recordless, update(1), &
135 & 1, mtc, 1, 1, 1, ntc(ng), 1, &
136 & tides(ng) % Tperiod)
137 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
138 END IF
139 END IF
140#endif
141
142#ifdef SSH_TIDES
143!
144! Tidal elevation amplitude and phase. In order to read data as a
145! function of tidal period, we need to reset the model time variables
146! temporarily.
147!
148 IF (lprocesstides(ng)) THEN
149 IF (iic(ng).eq.0) THEN
150 time_save=time(ng)
151 time(ng)=8640000.0_r8
152 tdays(ng)=time(ng)*sec2day
153!
154 CALL get_2dfld (ng, inlm, idtzam, tide(ng)%ncid, &
155# if defined PIO_LIB && defined DISTRIBUTE
156 & tide(ng)%pioFile, &
157# endif
158 & 1, tide(ng), update(1), &
159 & lbi, ubi, lbj, ubj, mtc, ntc(ng), &
160# ifdef MASKING
161 & grid(ng) % rmask, &
162# endif
163 & tides(ng) % SSH_Tamp)
164 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
165!
166 CALL get_2dfld (ng, inlm, idtzph, tide(ng)%ncid, &
167# if defined PIO_LIB && defined DISTRIBUTE
168 & tide(ng)%pioFile, &
169# endif
170 & 1, tide(ng), update(1), &
171 & lbi, ubi, lbj, ubj, mtc, ntc(ng), &
172# ifdef MASKING
173 & grid(ng) % rmask, &
174# endif
175 & tides(ng) % SSH_Tphase)
176 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
177!
178 time(ng)=time_save
179 tdays(ng)=time(ng)*sec2day
180 END IF
181 END IF
182#endif
183
184#ifdef UV_TIDES
185!
186! Tidal currents angle, phase, major and minor ellipse axis.
187!
188 IF (lprocesstides(ng)) THEN
189 IF (iic(ng).eq.0) THEN
190 time_save=time(ng)
191 time(ng)=8640000.0_r8
192 tdays(ng)=time(ng)*sec2day
193!
194 CALL get_2dfld (ng, inlm, idtvan, tide(ng)%ncid, &
195# if defined PIO_LIB && defined DISTRIBUTE
196 & tide(ng)%pioFile, &
197# endif
198 & 1, tide(ng), update(1), &
199 & lbi, ubi, lbj, ubj, mtc, ntc(ng), &
200# ifdef MASKING
201 & grid(ng) % rmask, &
202# endif
203 & tides(ng) % UV_Tangle)
204 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
205!
206 CALL get_2dfld (ng, inlm, idtvph, tide(ng)%ncid, &
207# if defined PIO_LIB && defined DISTRIBUTE
208 & tide(ng)%pioFile, &
209# endif
210 & 1, tide(ng), update(1), &
211 & lbi, ubi, lbj, ubj, mtc, ntc(ng), &
212# ifdef MASKING
213 & grid(ng) % rmask, &
214# endif
215 & tides(ng) % UV_Tphase)
216 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
217!
218 CALL get_2dfld (ng, inlm, idtvma, tide(ng)%ncid, &
219# if defined PIO_LIB && defined DISTRIBUTE
220 & tide(ng)%pioFile, &
221# endif
222 & 1, tide(ng), update(1), &
223 & lbi, ubi, lbj, ubj, mtc, ntc(ng), &
224# ifdef MASKING
225 & grid(ng) % rmask, &
226# endif
227 & tides(ng) % UV_Tmajor)
228 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
229!
230 CALL get_2dfld (ng, inlm, idtvmi, tide(ng)%ncid, &
231# if defined PIO_LIB && defined DISTRIBUTE
232 & tide(ng)%pioFile, &
233# endif
234 & 1, tide(ng), update(1), &
235 & lbi, ubi, lbj, ubj, mtc, ntc(ng), &
236# ifdef MASKING
237 & grid(ng) % rmask, &
238# endif
239 & tides(ng) % UV_Tminor)
240 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
241!
242 time(ng)=time_save
243 tdays(ng)=time(ng)*sec2day
244 END IF
245 END IF
246#endif
247
248#if defined AVERAGES && defined AVERAGES_DETIDE && \
249 (defined ssh_tides || defined uv_tides)
250!
251!-----------------------------------------------------------------------
252! If detiding and applicable, define additional variable to store
253! time-accumulated tide harmonics variables. This variable are
254! defined and written into input tide forcing NetCDF file.
255!-----------------------------------------------------------------------
256!
257 CALL def_tides (ng, ldeftide(ng))
258 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
259!
260!-----------------------------------------------------------------------
261! If restarting, read in time-accumulated tide harmonics variables.
262!-----------------------------------------------------------------------
263!
264 IF (.not.ldeftide(ng).and.(nrrec(ng).ne.0)) THEN
265!
266! For consistency, check time of written accumulate harmonics and
267! compare to current time.
268!
269 SELECT CASE (har(ng)%IOtype)
270 CASE (io_nf90)
271 CALL netcdf_get_time (ng, inlm, har(ng)%name, &
272 & vname(1,idtime), &
273 & rclock%DateNumber, htime, &
274 & ncid = har(ng)%ncid)
275
276# if defined PIO_LIB && defined DISTRIBUTE
277 CASE (io_pio)
278 CALL pio_netcdf_get_time (ng, inlm, har(ng)%name, &
279 & vname(1,idtime), &
280 & rclock%DateNumber, htime, &
281 & piofile = har(ng)%pioFile)
282# endif
283 END SELECT
284 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
285!
286 IF (time(ng).ne.htime) THEN
287 IF (master) THEN
288 WRITE (stdout,20) tdays(ng), htime*sec2day
289 END IF
290 exit_flag=2
291 ioerror=0
292 RETURN
293 END IF
294!
295! Number of time-acummulate tide harmonics.
296!
297 SELECT CASE (har(ng)%IOtype)
298 CASE (io_nf90)
299 CALL netcdf_get_ivar (ng, inlm, har(ng)%name, &
300 & 'Hcount', hcount(ng), &
301 & ncid = har(ng)%ncid)
302
303# if defined PIO_LIB && defined DISTRIBUT
304 CASE (io_pio)
305 CALL pio_netcdf_get_ivar (ng, inlm, har(ng)%name, &
306 & 'Hcount', hcount(ng), &
307 & piofile = har(ng)%pioFile)
308# endif
309 END SELECT
310 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
311!
312! Time-accumulated COS(omega(k)*t) harmonics.
313!
314 CALL get_ngfld (ng, inlm, idcosw, har(ng)%ncid, &
315# if defined PIO_LIB && defined DISTRIBUTE
316 & har(ng)%pioFile, &
317# endif
318 & 1, har(ng), recordless, update(1), &
319 & 1, mtc, 1, 1, 1, ntc(ng), 1, &
320 & tides(ng) % CosW_sum)
321 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
322!
323! Time-accumulated SIN(omega(k)*t) harmonics.
324!
325 CALL get_ngfld (ng, inlm, idsinw, har(ng)%ncid, &
326# if defined PIO_LIB && defined DISTRIBUTE
327 & har(ng)%pioFile, &
328# endif
329 & 1, har(ng), recordless, update(1), &
330 & 1, mtc, 1, 1, 1, ntc(ng), 1, &
331 & tides(ng) % SinW_sum)
332 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
333!
334! Time-accumulated COS(omega(k)*t)*COS(omega(k)*t) harmonics.
335!
336 CALL get_ngfld (ng, inlm, idcos2, har(ng)%ncid, &
337# if defined PIO_LIB && defined DISTRIBUTE
338 & har(ng)%pioFile, &
339# endif
340 & 1, har(ng), recordless, update(1), &
341 & 1, mtc, mtc, 1, 1, ntc(ng), ntc(ng), &
342 & tides(ng) % CosWCosW)
343 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
344!
345! Time-accumulated SIN(omega(k)*t)*SIN(omega(k)*t) harmonics.
346!
347 CALL get_ngfld (ng, inlm, idsin2, har(ng)%ncid, &
348# if defined PIO_LIB && defined DISTRIBUTE
349 & har(ng)%pioFile, &
350# endif
351 & 1, har(ng), recordless, update(1), &
352 & 1, mtc, mtc, 1, 1, ntc(ng), ntc(ng), &
353 & tides(ng) % SinWSinW)
354 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
355!
356! Time-accumulated SIN(omega(k)*t)*COS(omega(k)*t) harmonics.
357!
358 CALL get_ngfld (ng, inlm, idswcw, har(ng)%ncid, &
359# if defined PIO_LIB && defined DISTRIBUTE
360 & har(ng)%pioFile, &
361# endif
362 & 1, har(ng), recordless, update(1), &
363 & 1, mtc, mtc, 1, 1, ntc(ng), ntc(ng), &
364 & tides(ng) % SinWCosW)
365 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
366!
367! Time-accumulated free-surface tide harmonics.
368!
369 IF (aout(idfsud,ng)) THEN
370 SELECT CASE (har(ng)%IOtype)
371 CASE (io_nf90)
372 gtype=r3dvar
373 status=nf_fread3d(ng, inlm, har(ng)%name, &
374 & har(ng)%ncid, &
375 & vname(1,idfsuh), &
376 & har(ng)%Vid(idfsuh), &
377 & 0, gtype, vsize, &
378 & lbi, ubi, lbj, ubj, 0, 2*ntc(ng), &
379 & fscl, fmin, fmax, &
380# ifdef MASKING
381 & grid(ng) % rmask, &
382# endif
383 & tides(ng) % zeta_tide)
384
385# if defined PIO_LIB && defined DISTRIBUTE
386 CASE (io_pio)
387
388 piovar%vd => har(ng)%pioVar(idfsuh)%vd
389 IF (kind(tides(ng)%zeta_tide).eq.8) THEN
390 piovar%dkind=pio_double
391 iodesc => iodesc_dp_r2dhar(ng)
392 ELSE
393 piovar%dkind=pio_real
394 iodesc => iodesc_sp_r2dhar(ng)
395 END IF
396 piovar%gtype=r2dvar
397!
398 status=nf_fread3d(ng, inlm, har(ng)%name, &
399 & har(ng)%pioFile, &
400 & vname(1,idfsuh), &
401 & piovar, &
402 & 0, iodesc, vsize, &
403 & lbi, ubi, lbj, ubj, 0, 2*ntc(ng), &
404 & fscl, fmin, fmax, &
405# ifdef MASKING
406 & grid(ng) % rmask, &
407# endif
408 & tides(ng) % zeta_tide)
409
410# endif
411 END SELECT
412 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
413 IF (master) THEN
414 WRITE (stdout,10) trim(vname(1,idfsuh)), &
415 & trim(har(ng)%name)
416 END IF
417 exit_flag=2
418 ioerror=status
419 RETURN
420 ELSE
421 IF (master) THEN
422 WRITE (stdout,30) trim(vname(2,idfsuh)), fmin, fmax
423 END IF
424 END IF
425 END IF
426!
427! Time-accumulated 2D u-momentum tide harmonics.
428!
429 IF (aout(idu2dd,ng)) THEN
430 SELECT CASE (har(ng)%IOtype)
431 CASE (io_nf90)
432 gtype=u3dvar
433 status=nf_fread3d(ng, inlm, har(ng)%name, &
434 & har(ng)%ncid, &
435 & vname(1,idu2dh), &
436 & har(ng)%Vid(idu2dh), &
437 & 0, gtype, vsize, &
438 & lbi, ubi, lbj, ubj, 0, 2*ntc(ng), &
439 & fscl, fmin, fmax, &
440# ifdef MASKING
441 & grid(ng) % umask, &
442# endif
443 & tides(ng) % ubar_tide)
444
445# if defined PIO_LIB && defined DISTRIBUTE
446 CASE (io_pio)
447
448 piovar%vd => har(ng)%pioVar(idu2dh)%vd
449 IF (kind(tides(ng)%ubar_tide).eq.8) THEN
450 piovar%dkind=pio_double
451 iodesc => iodesc_dp_u2dhar(ng)
452 ELSE
453 piovar%dkind=pio_real
454 iodesc => iodesc_sp_u2dhar(ng)
455 END IF
456 piovar%gtype=u2dvar
457!
458 status=nf_fread3d(ng, inlm, har(ng)%name, &
459 & har(ng)%pioFile, &
460 & vname(1,idu2dh), &
461 & piovar, &
462 & 0, iodesc, vsize, &
463 & lbi, ubi, lbj, ubj, 0, 2*ntc(ng), &
464 & fscl, fmin, fmax, &
465# ifdef MASKING
466 & grid(ng) % umask, &
467# endif
468 & tides(ng) % ubar_tide)
469# endif
470 END SELECT
471 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
472 IF (master) THEN
473 WRITE (stdout,10) trim(vname(1,idu2dh)), &
474 & trim(har(ng)%name)
475 END IF
476 exit_flag=2
477 ioerror=status
478 RETURN
479 ELSE
480 IF (master) THEN
481 WRITE (stdout,30) trim(vname(2,idu2dh)), fmin, fmax
482 END IF
483 END IF
484 END IF
485!
486! Time-accumulated 2D v-momentum tide harmonics.
487!
488 IF (aout(idv2dd,ng)) THEN
489 SELECT CASE (har(ng)%IOtype)
490 CASE (io_nf90)
491 gtype=v3dvar
492 status=nf_fread3d(ng, inlm, har(ng)%name, &
493 & har(ng)%ncid, &
494 & vname(1,idv2dh), &
495 & har(ng)%Vid(idv2dh), &
496 & 0, gtype, vsize, &
497 & lbi, ubi, lbj, ubj, 0, 2*ntc(ng), &
498 & fscl, fmin, fmax, &
499# ifdef MASKING
500 & grid(ng) % vmask, &
501# endif
502 & tides(ng) % vbar_tide)
503
504# if defined PIO_LIB && defined DISTRIBUTE
505 CASE (io_pio)
506
507 piovar%vd => har(ng)%pioVar(idv2dh)%vd
508 IF (kind(tides(ng)%vbar_tide).eq.8) THEN
509 piovar%dkind=pio_double
510 iodesc => iodesc_dp_v2dhar(ng)
511 ELSE
512 piovar%dkind=pio_real
513 iodesc => iodesc_sp_v2dhar(ng)
514 END IF
515 piovar%gtype=v2dvar
516!
517 status=nf_fread3d(ng, inlm, har(ng)%name, &
518 & har(ng)%pioFile, &
519 & vname(1,idv2dh), &
520 & piovar, &
521 & 0, iodesc, vsize, &
522 & lbi, ubi, lbj, ubj, 0, 2*ntc(ng), &
523 & fscl, fmin, fmax, &
524# ifdef MASKING
525 & grid(ng) % vmask, &
526# endif
527 & tides(ng) % vbar_tide)
528# endif
529 END SELECT
530 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
531 IF (master) THEN
532 WRITE (stdout,10) trim(vname(1,idv2dh)), &
533 & trim(har(ng)%name)
534 END IF
535 exit_flag=2
536 ioerror=status
537 RETURN
538 ELSE
539 IF (master) THEN
540 WRITE (stdout,30) trim(vname(2,idv2dh)), fmin, fmax
541 END IF
542 END IF
543 END IF
544
545# ifdef SOLVE3D
546!
547! Time-accumulated 3D u-momentum tide harmonics.
548!
549 IF (aout(idu3dd,ng)) THEN
550 SELECT CASE (har(ng)%IOtype)
551 CASE (io_nf90)
552 gtype=u3dvar
553 status=nf_fread4d(ng, inlm, har(ng)%name, &
554 & har(ng)%ncid, &
555 & vname(1,idu3dh), &
556 & har(ng)%Vid(idu3dh), &
557 & 0, gtype, vsize, &
558 & lbi, ubi, lbj, ubj, 1, n(ng), &
559 & 0, 2*ntc(ng), &
560 & fscl, fmin, fmax, &
561# ifdef MASKING
562 & grid(ng) % umask, &
563# endif
564 & tides(ng) % u_tide)
565
566# if defined PIO_LIB && defined DISTRIBUTE
567 CASE (io_pio)
568
569 piovar%vd => har(ng)%pioVar(idu3dh)%vd
570 IF (kind(tides(ng)%u_tide).eq.8) THEN
571 piovar%dkind=pio_double
572 iodesc => iodesc_dp_u3dhar(ng)
573 ELSE
574 piovar%dkind=pio_real
575 iodesc => iodesc_sp_u3dhar(ng)
576 END IF
577 piovar%gtype=u3dvar
578!
579 status=nf_fread4d(ng, inlm, har(ng)%name, &
580 & har(ng)%pioFile, &
581 & vname(1,idu3dh), &
582 & piovar, &
583 & 0, iodesc, vsize, &
584 & lbi, ubi, lbj, ubj, 1, n(ng), &
585 & 0, 2*ntc(ng), &
586 & fscl, fmin, fmax, &
587# ifdef MASKING
588 & grid(ng) % umask, &
589# endif
590 & tides(ng) % u_tide)
591# endif
592 END SELECT
593 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
594 IF (master) THEN
595 WRITE (stdout,10) trim(vname(1,idu3dh)), &
596 & trim(har(ng)%name)
597 END IF
598 exit_flag=2
599 ioerror=status
600 RETURN
601 ELSE
602 IF (master) THEN
603 WRITE (stdout,30) trim(vname(2,idu3dh)), fmin, fmax
604 END IF
605 END IF
606 END IF
607!
608! Time-accumulated 3D v-momentum tide harmonics.
609!
610 IF (aout(idv3dd,ng)) THEN
611 SELECT CASE (har(ng)%IOtype)
612 CASE (io_nf90)
613 gtype=v3dvar
614 status=nf_fread4d(ng, inlm, har(ng)%name, &
615 & har(ng)%ncid, &
616 & vname(1,idv3dh), &
617 & har(ng)%Vid(idv3dh), &
618 & 0, gtype, vsize, &
619 & lbi, ubi, lbj, ubj, 1, n(ng), &
620 & 0, 2*ntc(ng), &
621 & fscl, fmin, fmax, &
622# ifdef MASKING
623 & grid(ng) % vmask, &
624# endif
625 & tides(ng) % v_tide)
626
627# if defined PIO_LIB && defined DISTRIBUTE
628 CASE (io_pio)
629
630 piovar%vd => har(ng)%pioVar(idv3dh)%vd
631 IF (kind(tides(ng)%v_tide).eq.8) THEN
632 piovar%dkind=pio_double
633 iodesc => iodesc_dp_v3dhar(ng)
634 ELSE
635 piovar%dkind=pio_real
636 iodesc => iodesc_sp_v3dhar(ng)
637 END IF
638 piovar%gtype=v3dvar
639!
640 status=nf_fread4d(ng, inlm, har(ng)%name, &
641 & har(ng)%pioFile, &
642 & vname(1,idv3dh), &
643 & piovar, &
644 & 0, iodesc, vsize, &
645 & lbi, ubi, lbj, ubj, 1, n(ng), &
646 & 0, 2*ntc(ng), &
647 & fscl, fmin, fmax, &
648# ifdef MASKING
649 & grid(ng) % vmask, &
650# endif
651 & tides(ng) % v_tide)
652# endif
653 END SELECT
654 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
655 IF (master) THEN
656 WRITE (stdout,10) trim(vname(1,idv3dh)), &
657 & trim(har(ng)%name)
658 END IF
659 exit_flag=2
660 ioerror=status
661 RETURN
662 ELSE
663 IF (master) THEN
664 WRITE (stdout,30) trim(vname(2,idv3dh)), fmin, fmax
665 END IF
666 END IF
667 END IF
668!
669! Time-accumulated temperature and salinity tide harmonics.
670!
671 DO itrc=1,nat
672 IF (aout(idtrcd(itrc),ng)) THEN
673 SELECT CASE (har(ng)%IOtype)
674 CASE (io_nf90)
675 gtype=r3dvar
676 status=nf_fread4d(ng, inlm, har(ng)%name, &
677 & har(ng)%ncid, &
678 & vname(1,idtrch(itrc)), &
679 & har(ng)%Vid(idtrch(itrc)), &
680 & 0, gtype, vsize, &
681 & lbi, ubi, lbj, ubj, 1, n(ng), &
682 & 0, 2*ntc(ng), &
683 & fscl, fmin, fmax, &
684# ifdef MASKING
685 & grid(ng) % rmask, &
686# endif
687 & tides(ng) % t_tide(:,:,:,:,itrc))
688
689# if defined PIO_LIB && defined DISTRIBUTE
690 CASE (io_pio)
691
692 piovar%vd => har(ng)%pioVar(idtrch(itrc))%vd
693 IF (kind(tides(ng)%t_tide).eq.8) THEN
694 piovar%dkind=pio_double
695 iodesc => iodesc_dp_r3dhar(ng)
696 ELSE
697 piovar%dkind=pio_real
698 iodesc => iodesc_sp_r3dhar(ng)
699 END IF
700 piovar%gtype=r3dvar
701!
702 status=nf_fread4d(ng, inlm, har(ng)%name, &
703 & har(ng)%pioFile, &
704 & vname(1,idtrch(itrc)), &
705 & piovar, &
706 & 0, iodesc, vsize, &
707 & lbi, ubi, lbj, ubj, 1, n(ng), &
708 & 0, 2*ntc(ng), &
709 & fscl, fmin, fmax, &
710# ifdef MASKING
711 & grid(ng) % rmask, &
712# endif
713 & tides(ng) % t_tide(:,:,:,:,itrc))
714# endif
715 END SELECT
716 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
717 IF (master) THEN
718 WRITE (stdout,10) trim(vname(1,idtrch(itrc))), &
719 & trim(har(ng)%name)
720 END IF
721 exit_flag=2
722 ioerror=status
723 RETURN
724 ELSE
725 IF (master) THEN
726 WRITE (stdout,30) trim(vname(2,idtrch(itrc))), &
727 & fmin, fmax
728 END IF
729 END IF
730 END IF
731 END DO
732# endif
733 END IF
734#endif
735
736#ifndef ANA_PSOURCE
737!
738!-----------------------------------------------------------------------
739! Read in point Sources/Sinks position, direction, special flag, and
740! mass transport nondimensional shape profile. Point sources are at
741! U- and V-points.
742!-----------------------------------------------------------------------
743!
744 IF ((iic(ng).eq.0).and. &
745 & (luvsrc(ng).or.lwsrc(ng).or.any(ltracersrc(:,ng)))) THEN
746 CALL get_ngfld (ng, inlm, idrxpo, ssf(ng)%ncid, &
747# if defined PIO_LIB && defined DISTRIBUTE
748 & ssf(ng)%pioFile, &
749# endif
750 & 1, ssf(ng), recordless, update(1), &
751 & 1, nsrc(ng), 1, 1, 1, nsrc(ng), 1, &
752 & sources(ng) % Xsrc)
753 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
754
755 CALL get_ngfld (ng, inlm, idrepo, ssf(ng)%ncid, &
756# if defined PIO_LIB && defined DISTRIBUTE
757 & ssf(ng)%pioFile, &
758# endif
759 & 1, ssf(ng), recordless, update(1), &
760 & 1, nsrc(ng), 1, 1, 1, nsrc(ng), 1, &
761 & sources(ng) % Ysrc)
762 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
763
764 CALL get_ngfld (ng, inlm, idrdir, ssf(ng)%ncid, &
765# if defined PIO_LIB && defined DISTRIBUTE
766 & ssf(ng)%pioFile, &
767# endif
768 & 1, ssf(ng), recordless, update(1), &
769 & 1, nsrc(ng), 1, 1, 1, nsrc(ng), 1, &
770 & sources(ng) % Dsrc)
771 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
772
773# ifdef SOLVE3D
774 CALL get_ngfld (ng, inlm, idrvsh, ssf(ng)%ncid, &
775# if defined PIO_LIB && defined DISTRIBUTE
776 & ssf(ng)%pioFile, &
777# endif
778 & 1, ssf(ng), recordless, update(1), &
779 & 1, nsrc(ng), n(ng), 1, 1, nsrc(ng), n(ng), &
780 & sources(ng) % Qshape)
781 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
782# endif
783
784 DO is=1,nsrc(ng)
785 sources(ng)%Isrc(is)= &
786 & max(1,min(nint(sources(ng)%Xsrc(is)),lm(ng)+1))
787 sources(ng)%Jsrc(is)= &
788 & max(1,min(nint(sources(ng)%Ysrc(is)),mm(ng)+1))
789 END DO
790 END IF
791#endif
792
793#ifdef RED_TIDE
794!
795!-----------------------------------------------------------------------
796! Red tides.
797!-----------------------------------------------------------------------
798!
799! Read initial dynoflagellate bottom cyst concentration.
800!
801 CALL get_2dfld (ng, inlm, idcyst, frcncid(idcyst,ng), &
802# if defined PIO_LIB && defined DISTRIBUTE
803 & frcpiofile(idcyst,ng), &
804# endif
805 & nffiles(ng), frc(1,ng), update(1), &
806 & lbi, ubi, lbj, ubj, 1, 1, &
807# ifdef MASKING
808 & grid(ng) % rmask, &
809# endif
810 & ocean(ng) % CystIni)
811 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
812 IF (.not.linfo(1,idcyst,ng)) THEN ! if not gridded,
813 ocean(ng) % CystIni = fpoint(1,idcyst,ng) ! load point data
814 END IF
815#endif
816
817#ifdef PROFILE
818!
819!-----------------------------------------------------------------------
820! Turn off input data time wall clock.
821!-----------------------------------------------------------------------
822!
823 CALL wclock_off (ng, inlm, 3, __line__, myfile)
824#endif
825
826#if defined AVERAGES && defined AVERAGES_DETIDE && \
827 (defined ssh_tides || defined uv_tides)
828!
829 10 FORMAT (/,' GET_IDATA - error while reading variable: ',a, &
830 & /,13x,'in input NetCDF file: ',a)
831 20 FORMAT (/,' GET_IDATA - incosistent restart and harmonics time:', &
832 & /,13x,f15.4,2x,f15.4)
833 30 FORMAT (16x,'- ',a,/,19x,'(Min = ',1p,e15.8, &
834 & ' Max = ',1p,e15.8,')')
835#endif
836!
837 RETURN
subroutine get_2dfld(ng, model, ifield, ncid, piofile, nfiles, s, update, lbi, ubi, lbj, ubj, iout, irec, fmask, fout)
Definition get_2dfld.F:12
subroutine get_ngfld(ng, model, ifield, ncid, piofile, nfiles, s, recordless, update, lbi, ubi, ubj, ubk, istr, iend, jrec, fout)
Definition get_ngfld.F:9
subroutine, public def_tides(ng, ldef)
Definition def_tides.F:45
integer idcyst
type(t_grid), dimension(:), allocatable grid
Definition mod_grid.F:365
integer ioerror
type(t_io), dimension(:), allocatable ssf
type(t_io), dimension(:), allocatable har
type(t_io), dimension(:,:), allocatable frc
type(t_io), dimension(:), allocatable tide
integer stdout
character(len=256) sourcefile
integer, dimension(:), allocatable nffiles
integer, parameter io_nf90
Definition mod_ncparam.F:95
integer idu2dh
type(file_desc_t), dimension(:,:), pointer frcpiofile
integer idu3dd
integer idv3dd
integer idsin2
integer idtzph
integer idtper
logical, dimension(:,:,:), allocatable linfo
integer, parameter io_pio
Definition mod_ncparam.F:96
integer idu3dh
integer idv3dh
integer idcos2
integer, dimension(:,:), allocatable frcncid
integer idfsud
integer idtvmi
integer idrepo
integer idtzam
integer idtvph
real(dp), dimension(:,:,:), allocatable fpoint
integer idswcw
integer, dimension(:), allocatable idtrcd
integer, dimension(:), allocatable idtrch
integer idtvan
integer idcosw
integer idrvsh
integer idv2dd
character(len=maxlen), dimension(6, 0:nv) vname
integer idtime
integer idrxpo
integer idtvma
integer idv2dh
integer idu2dd
integer idfsuh
integer idsinw
logical, dimension(:,:), allocatable aout
integer idrdir
type(t_ocean), dimension(:), allocatable ocean
Definition mod_ocean.F:351
logical master
integer nat
Definition mod_param.F:499
integer, parameter inlm
Definition mod_param.F:662
integer, dimension(:), allocatable n
Definition mod_param.F:479
integer mtc
Definition mod_param.F:564
integer, parameter r3dvar
Definition mod_param.F:721
integer, parameter u3dvar
Definition mod_param.F:722
integer, dimension(:), allocatable lm
Definition mod_param.F:455
integer, parameter u2dvar
Definition mod_param.F:718
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
type(io_desc_t), dimension(:), pointer iodesc_dp_r2dhar
type(io_desc_t), dimension(:), pointer iodesc_dp_u3dhar
type(io_desc_t), dimension(:), pointer iodesc_dp_v3dhar
type(io_desc_t), dimension(:), pointer iodesc_sp_v2dhar
type(io_desc_t), dimension(:), pointer iodesc_dp_v2dhar
type(io_desc_t), dimension(:), pointer iodesc_sp_r3dhar
type(io_desc_t), dimension(:), pointer iodesc_sp_r2dhar
type(io_desc_t), dimension(:), pointer iodesc_sp_u2dhar
type(io_desc_t), dimension(:), pointer iodesc_sp_v3dhar
type(io_desc_t), dimension(:), pointer iodesc_sp_u3dhar
type(io_desc_t), dimension(:), pointer iodesc_dp_u2dhar
type(io_desc_t), dimension(:), pointer iodesc_dp_r3dhar
logical, dimension(:), allocatable luvsrc
integer, dimension(:), allocatable nrrec
logical, dimension(:,:), allocatable ltracersrc
integer, dimension(:), allocatable iic
logical, dimension(:), allocatable lprocesstides
logical, dimension(:), allocatable ldeftide
real(dp), dimension(:), allocatable tdays
type(t_clock) rclock
logical, dimension(:), allocatable lwsrc
real(dp), parameter sec2day
integer, dimension(:), allocatable hcount
integer exit_flag
real(dp), dimension(:), allocatable time
integer noerror
type(t_sources), dimension(:), allocatable sources
Definition mod_sources.F:90
integer, dimension(:), allocatable nsrc
Definition mod_sources.F:97
integer, dimension(:), allocatable ntc
type(t_tides), dimension(:), allocatable tides
Definition mod_tides.F:133
logical function, public founderror(flag, noerr, line, routine)
Definition strings.F:52
recursive subroutine wclock_off(ng, model, region, line, routine)
Definition timers.F:148
recursive subroutine wclock_on(ng, model, region, line, routine)
Definition timers.F:3

References mod_ncparam::aout, def_tides_mod::def_tides(), mod_scalars::exit_flag, strings_mod::founderror(), mod_ncparam::fpoint, mod_iounits::frc, mod_ncparam::frcncid, mod_ncparam::frcpiofile, get_2dfld(), get_ngfld(), mod_grid::grid, mod_iounits::har, mod_scalars::hcount, mod_ncparam::idcos2, mod_ncparam::idcosw, mod_biology::idcyst, mod_ncparam::idfsud, mod_ncparam::idfsuh, mod_ncparam::idrdir, mod_ncparam::idrepo, mod_ncparam::idrvsh, mod_ncparam::idrxpo, mod_ncparam::idsin2, mod_ncparam::idsinw, mod_ncparam::idswcw, mod_ncparam::idtime, mod_ncparam::idtper, mod_ncparam::idtrcd, mod_ncparam::idtrch, mod_ncparam::idtvan, mod_ncparam::idtvma, mod_ncparam::idtvmi, mod_ncparam::idtvph, mod_ncparam::idtzam, mod_ncparam::idtzph, mod_ncparam::idu2dd, mod_ncparam::idu2dh, mod_ncparam::idu3dd, mod_ncparam::idu3dh, mod_ncparam::idv2dd, mod_ncparam::idv2dh, mod_ncparam::idv3dd, mod_ncparam::idv3dh, mod_scalars::iic, mod_param::inlm, mod_ncparam::io_nf90, mod_ncparam::io_pio, mod_pio_netcdf::iodesc_dp_r2dhar, mod_pio_netcdf::iodesc_dp_r3dhar, mod_pio_netcdf::iodesc_dp_u2dhar, mod_pio_netcdf::iodesc_dp_u3dhar, mod_pio_netcdf::iodesc_dp_v2dhar, mod_pio_netcdf::iodesc_dp_v3dhar, mod_pio_netcdf::iodesc_sp_r2dhar, mod_pio_netcdf::iodesc_sp_r3dhar, mod_pio_netcdf::iodesc_sp_u2dhar, mod_pio_netcdf::iodesc_sp_u3dhar, mod_pio_netcdf::iodesc_sp_v2dhar, mod_pio_netcdf::iodesc_sp_v3dhar, mod_iounits::ioerror, mod_scalars::ldeftide, mod_ncparam::linfo, mod_param::lm, mod_scalars::lprocesstides, mod_scalars::ltracersrc, mod_scalars::luvsrc, mod_scalars::lwsrc, mod_parallel::master, mod_param::mm, mod_param::mtc, mod_param::n, mod_param::nat, mod_iounits::nffiles, mod_scalars::noerror, mod_scalars::nrrec, mod_sources::nsrc, mod_stepping::ntc, mod_ocean::ocean, mod_param::r2dvar, mod_param::r3dvar, mod_scalars::rclock, mod_scalars::sec2day, mod_iounits::sourcefile, mod_sources::sources, mod_iounits::ssf, mod_iounits::stdout, mod_scalars::tdays, mod_iounits::tide, mod_tides::tides, mod_scalars::time, mod_param::u2dvar, mod_param::u3dvar, mod_param::v2dvar, mod_param::v3dvar, mod_ncparam::vname, wclock_off(), and wclock_on().

Referenced by initial(), and roms_kernel_mod::nlm_initial().

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