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

Go to the source code of this file.

Functions/Subroutines

subroutine tl_get_data (ng)
 

Function/Subroutine Documentation

◆ tl_get_data()

subroutine tl_get_data ( integer, intent(in) ng)

Definition at line 3 of file tl_get_data.F.

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 routine reads in forcing, climatology and other data from !
13! NetCDF files. If there is more than one time-record, data is !
14! loaded into global two-time record arrays. The interpolation !
15! is carried elsewhere. !
16! !
17! Currently, this routine is only executed in serial mode by the !
18! main thread. !
19! !
20!=======================================================================
21!
22 USE mod_param
23 USE mod_boundary
24# ifdef FORWARD_READ
25 USE mod_coupling
26# endif
27 USE mod_clima
28 USE mod_forces
29 USE mod_grid
30 USE mod_iounits
31 USE mod_mixing
32 USE mod_ncparam
33# ifdef FORWARD_READ
34 USE mod_ocean
35# endif
36 USE mod_scalars
37 USE mod_sources
38 USE mod_stepping
39!
40 USE strings_mod, ONLY : founderror
41!
42 implicit none
43!
44! Imported variable declarations.
45!
46 integer, intent(in) :: ng
47!
48! Local variable declarations.
49!
50 logical, save :: recordless = .false.
51
52 logical, dimension(3) :: update = &
53 & (/ .FALSE., .FALSE., .FALSE. /)
54!
55 integer :: ILB, IUB, JLB, JUB
56 integer :: LBi, UBi, LBj, UBj
57 integer :: i, ic, my_tile
58
59# ifdef FORWARD_MIXING
60!
61 real(r8) :: scale
62# endif
63!
64 character (len=*), parameter :: MyFile = &
65 & __FILE__
66!
67! Lower and upper bounds for nontiled (global values) boundary arrays.
68!
69 my_tile=-1 ! for global values
70 ilb=bounds(ng)%LBi(my_tile)
71 iub=bounds(ng)%UBi(my_tile)
72 jlb=bounds(ng)%LBj(my_tile)
73 jub=bounds(ng)%UBj(my_tile)
74!
75! Lower and upper bounds for tiled arrays.
76!
77 lbi=lbound(grid(ng)%h,dim=1)
78 ubi=ubound(grid(ng)%h,dim=1)
79 lbj=lbound(grid(ng)%h,dim=2)
80 ubj=ubound(grid(ng)%h,dim=2)
81
82# ifdef PROFILE
83!
84!-----------------------------------------------------------------------
85! Turn on input data time wall clock.
86!-----------------------------------------------------------------------
87!
88 CALL wclock_on (ng, itlm, 3, __line__, myfile)
89# endif
90!
91!=======================================================================
92! Read in forcing data from FORCING NetCDF file.
93!=======================================================================
94
95# ifndef ANA_PSOURCE
96!
97!-----------------------------------------------------------------------
98! Point Sources/Sinks time dependent data.
99!-----------------------------------------------------------------------
100!
101! Point Source/Sink vertically integrated mass transport.
102!
103 IF (luvsrc(ng).or.lwsrc(ng)) THEN
104 CALL get_ngfld (ng, itlm, idrtra, ssf(ng)%ncid, &
105# if defined PIO_LIB && defined DISTRIBUTE
106 & ssf(ng)%pioFile, &
107# endif
108 & 1, ssf(ng), recordless, update(1), &
109 & 1, nsrc(ng), 1, 2, 1, nsrc(ng), 1, &
110 & sources(ng) % QbarG)
111 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
112 END IF
113
114# ifdef SOLVE3D
115!
116! Tracer Sources/Sinks.
117!
118 DO i=1,nt(ng)
119 IF (ltracersrc(i,ng)) THEN
120 CALL get_ngfld (ng, itlm, idrtrc(i), ssf(ng)%ncid, &
121# if defined PIO_LIB && defined DISTRIBUTE
122 & ssf(ng)%pioFile, &
123# endif
124 & 1, ssf(ng), recordless, update(1), &
125 & 1, nsrc(ng), n(ng), 2, 1, nsrc(ng), n(ng), &
126 & sources(ng) % TsrcG(:,:,:,i))
127 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
128 END IF
129 END DO
130# endif
131# endif
132
133# if !defined ANA_WINDS && \
134 ((defined bulk_fluxes && !defined FORWARD_FLUXES) || \
135 defined ecosim)
136!
137!-----------------------------------------------------------------------
138! Surface wind components.
139!-----------------------------------------------------------------------
140!
141 CALL get_2dfld (ng, itlm, iduair, frcncid(iduair,ng), &
142# if defined PIO_LIB && defined DISTRIBUTE
143 & frcpiofile(iduair,ng), &
144# endif
145 & nffiles(ng), frc(1,ng), update(1), &
146 & lbi, ubi, lbj, ubj, 2, 1, &
147# ifdef MASKING
148 & grid(ng) % rmask, &
149# endif
150 & forces(ng) % UwindG)
151 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
152!
153 CALL get_2dfld (ng , itlm, idvair, frcncid(idvair,ng), &
154# if defined PIO_LIB && defined DISTRIBUTE
155 & frcpiofile(idvair,ng), &
156# endif
157 & nffiles(ng), frc(1,ng), update(1), &
158 & lbi, ubi, lbj, ubj, 2, 1, &
159# ifdef MASKING
160 & grid(ng) % rmask, &
161# endif
162 & forces(ng) % VwindG)
163 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
164# endif
165
166# ifndef FRC_COUPLING
167# if !defined ANA_SMFLUX && \
168 !defined BULK_FLUXES && !defined FORWARD_FLUXES
169!
170!-----------------------------------------------------------------------
171! Surface wind stress components from input FRC file.
172!-----------------------------------------------------------------------
173!
174 CALL get_2dfld (ng, itlm, idusms, frcncid(idusms,ng), &
175# if defined PIO_LIB && defined DISTRIBUTE
176 & frcpiofile(idusms,ng), &
177# endif
178 & nffiles(ng), frc(1,ng), update(1), &
179 & lbi, ubi, lbj, ubj, 2, 1, &
180# ifdef MASKING
181 & grid(ng) % umask, &
182# endif
183 & forces(ng) % sustrG)
184 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
185!
186 CALL get_2dfld (ng, itlm, idvsms, frcncid(idvsms,ng), &
187# if defined PIO_LIB && defined DISTRIBUTE
188 & frcpiofile(idvsms,ng), &
189# endif
190 & nffiles(ng), frc(1,ng), update(1), &
191 & lbi, ubi, lbj, ubj, 2, 1, &
192# ifdef MASKING
193 & grid(ng) % vmask, &
194# endif
195 & forces(ng) % svstrG)
196 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
197# endif
198# endif
199
200# if defined FORWARD_FLUXES || defined FRC_COUPLING
201!
202!-----------------------------------------------------------------------
203! Surface wind stress components from NLM forward file.
204!-----------------------------------------------------------------------
205!
206 CALL get_2dfld (ng, itlm, idusms, blk(ng)%ncid, &
207# if defined PIO_LIB && defined DISTRIBUTE
208 & blk(ng)%pioFile, &
209# endif
210 & 1, blk(ng), update(1), &
211 & lbi, ubi, lbj, ubj, 2, 1, &
212# ifdef MASKING
213 & grid(ng) % umask, &
214# endif
215 & forces(ng) % sustrG)
216 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
217!
218 CALL get_2dfld (ng, itlm, idvsms, blk(ng)%ncid, &
219# if defined PIO_LIB && defined DISTRIBUTE
220 & blk(ng)%pioFile, &
221# endif
222 & 1, blk(ng), update(1), &
223 & lbi, ubi, lbj, ubj, 2, 1, &
224# ifdef MASKING
225 & grid(ng) % vmask, &
226# endif
227 & forces(ng) % svstrG)
228 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
229# endif
230
231# if !defined ANA_PAIR && \
232 (defined atm_press || defined bulk_fluxes || \
233 defined ecosim)
234# if !(defined FRC_COUPLING || defined FORWARD_FLUXES)
235!
236!-----------------------------------------------------------------------
237! Surface air pressure from input FRC file.
238!-----------------------------------------------------------------------
239!
240 CALL get_2dfld (ng, itlm, idpair, frcncid(idpair,ng), &
241# if defined PIO_LIB && defined DISTRIBUTE
242 & frcpiofile(idpair,ng), &
243# endif
244 & nffiles(ng), frc(1,ng), update(1), &
245 & lbi, ubi, lbj, ubj, 2, 1, &
246# ifdef MASKING
247 & grid(ng) % rmask, &
248# endif
249 & forces(ng) % PairG)
250 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
251
252# else
253!
254!-----------------------------------------------------------------------
255! Surface air pressure from NLM forward file.
256!-----------------------------------------------------------------------
257!
258 CALL get_2dfld (ng, itlm, idpair, blk(ng)%ncid, &
259# if defined PIO_LIB && defined DISTRIBUTE
260 & blk(ng)%pioFile, &
261# endif
262 & 1, blk(ng), update(1), &
263 & lbi, ubi, lbj, ubj, 2, 1, &
264# ifdef MASKING
265 & grid(ng) % rmask, &
266# endif
267 & forces(ng) % PairG)
268 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
269# endif
270# endif
271
272# if !defined ANA_WWAVE && defined WAVE_DATA
273!
274!-----------------------------------------------------------------------
275! Surface wind induced wave properties.
276!-----------------------------------------------------------------------
277!
278# ifdef WAVES_DIR
279 CALL get_2dfld (ng, itlm, idwdir, frcncid(idwdir,ng), &
280# if defined PIO_LIB && defined DISTRIBUTE
281 & frcpiofile(idwdir,ng), &
282# endif
283 & nffiles(ng), frc(1,ng), update(1), &
284 & lbi, ubi, lbj, ubj, 2, 1, &
285# ifdef MASKING
286 & grid(ng) % rmask, &
287# endif
288 & forces(ng) % DwaveG)
289 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
290# endif
291
292# ifdef WAVES_HEIGHT
293!
294 CALL get_2dfld (ng, itlm, idwamp, frcncid(idwamp,ng), &
295# if defined PIO_LIB && defined DISTRIBUTE
296 & frcpiofile(idwamp,ng), &
297# endif
298 & nffiles(ng), frc(1,ng), update(1), &
299 & lbi, ubi, lbj, ubj, 2, 1, &
300# ifdef MASKING
301 & grid(ng) % rmask, &
302# endif
303 & forces(ng) % HwaveG)
304 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
305# endif
306
307# ifdef WAVES_LENGTH
308!
309 CALL get_2dfld (ng, itlm, idwlen, frcncid(idwlen,ng), &
310# if defined PIO_LIB && defined DISTRIBUTE
311 & frcpiofile(idwlen,ng), &
312# endif
313 & nffiles(ng), frc(1,ng), update(1), &
314 & lbi, ubi, lbj, ubj, 2, 1, &
315# ifdef MASKING
316 & grid(ng) % rmask, &
317# endif
318 & forces(ng) % LwaveG)
319 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
320# endif
321
322# ifdef WAVES_TOP_PERIOD
323!
324 CALL get_2dfld (ng, itlm, idwptp, frcncid(idwptp,ng), &
325# if defined PIO_LIB && defined DISTRIBUTE
326 & frcpiofile(idwptp,ng), &
327# endif
328 & nffiles(ng), frc(1,ng), update(1), &
329 & lbi, ubi, lbj, ubj, 2, 1, &
330# ifdef MASKING
331 & grid(ng) % rmask, &
332# endif
333 & forces(ng) % Pwave_topG)
334 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
335# endif
336
337# ifdef WAVES_BOT_PERIOD
338!
339 CALL get_2dfld (ng, itlm, idwpbt, frcncid(idwpbt,ng), &
340# if defined PIO_LIB && defined DISTRIBUTE
341 & frcpiofile(idwpbt,ng), &
342# endif
343 & nffiles(ng), frc(1,ng), update(1), &
344 & lbi, ubi, lbj, ubj, 2, 1, &
345# ifdef MASKING
346 & grid(ng) % rmask, &
347# endif
348 & forces(ng) % Pwave_botG)
349 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
350# endif
351
352# if defined WAVES_UB
353!
354 CALL get_2dfld (ng, itlm, idworb, frcncid(idworb,ng), &
355# if defined PIO_LIB && defined DISTRIBUTE
356 & frcpiofile(idworb,ng), &
357# endif
358 & nffiles(ng), frc(1,ng), update(1), &
359 & lbi, ubi, lbj, ubj, 2, 1, &
360# ifdef MASKING
361 & grid(ng) % rmask, &
362# endif
363 & forces(ng) % Ub_swanG)
364 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
365# endif
366
367# if defined TKE_WAVEDISS
368!
369 CALL get_2dfld (ng, itlm, idwdis, frcncid(idwdis,ng), &
370# if defined PIO_LIB && defined DISTRIBUTE
371 & frcpiofile(idwdis,ng), &
372# endif
373 & nffiles(ng), frc(1,ng), update(1), &
374 & lbi, ubi, lbj, ubj, 2, 1, &
375# ifdef MASKING
376 & grid(ng) % rmask, &
377# endif
378 & forces(ng) % Wave_dissipG)
379 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
380# endif
381
382# if defined SVENDSEN_ROLLER
383!
384 CALL get_2dfld (ng, itlm, idwbrk, frcncid(idwbrk,ng), &
385# if defined PIO_LIB && defined DISTRIBUTE
386 & frcpiofile(idwbrk,ng), &
387# endif
388 & nffiles(ng), frc(1,ng), update(1), &
389 & lbi, ubi, lbj, ubj, 2, 1, &
390# ifdef MASKING
391 & grid(ng) % rmask, &
392# endif
393 & forces(ng) % Wave_breakG)
394 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
395# endif
396# endif
397
398# ifdef SOLVE3D
399
400# if !defined ANA_CLOUD && defined CLOUDS
401!
402!-----------------------------------------------------------------------
403! Cloud fraction.
404!-----------------------------------------------------------------------
405!
406 CALL get_2dfld (ng, itlm, idcfra, frcncid(idcfra,ng), &
407# if defined PIO_LIB && defined DISTRIBUTE
408 & frcpiofile(idcfra,ng), &
409# endif
410 & nffiles(ng), frc(1,ng), update(1), &
411 & lbi, ubi, lbj, ubj, 2, 1, &
412# ifdef MASKING
413 & grid(ng) % rmask, &
414# endif
415 & forces(ng) % cloudG)
416 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
417# endif
418
419# if defined SHORTWAVE && !defined ANA_SRFLUX
420# ifdef FORWARD_FLUXES
421!
422!-----------------------------------------------------------------------
423! Surface solar shortwave radiation from NLM forward file.
424!-----------------------------------------------------------------------
425!
426 CALL get_2dfld (ng, itlm, idsrad, blk(ng)%ncid, &
427# if defined PIO_LIB && defined DISTRIBUTE
428 & blk(ng)%pioFile, &
429# endif
430 & 1, blk(ng), update(1), &
431 & lbi, ubi, lbj, ubj, 2, 1, &
432# ifdef MASKING
433 & grid(ng) % rmask, &
434# endif
435 & forces(ng) % srflxG)
436 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
437
438# else
439
440# if !(defined BULK_FLUXES || defined FRC_COUPLING)
441!
442!-----------------------------------------------------------------------
443! Surface solar shortwave radiation.
444!-----------------------------------------------------------------------
445!
446 CALL get_2dfld (ng, itlm, idsrad, frcncid(idsrad,ng), &
447# if defined PIO_LIB && defined DISTRIBUTE
448 & frcpiofile(idsrad,ng), &
449# endif
450 & nffiles(ng), frc(1,ng), update(1), &
451 & lbi, ubi, lbj, ubj, 2, 1, &
452# ifdef MASKING
453 & grid(ng) % rmask, &
454# endif
455 & forces(ng) % srflxG)
456 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
457# endif
458
459# endif
460# endif
461
462# if (defined BULK_FLUXES && !defined FORWARD_FLUXES) && \
463 !defined LONGWAVE && !defined LONGWAVE_OUT
464!
465!-----------------------------------------------------------------------
466! Surface net longwave radiation.
467!-----------------------------------------------------------------------
468!
469 CALL get_2dfld (ng, itlm, idlrad, frcncid(idlrad,ng), &
470# if defined PIO_LIB && defined DISTRIBUTE
471 & frcpiofile(idlrad,ng), &
472# endif
473 & nffiles(ng), frc(1,ng), update(1), &
474 & lbi, ubi, lbj, ubj, 2, 1, &
475# ifdef MASKING
476 & grid(ng) % rmask, &
477# endif
478 & forces(ng) % lrflxG)
479 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
480# endif
481
482# if (defined BULK_FLUXES && !defined FORWARD_FLUXES) && \
483 defined longwave_out
484!
485!-----------------------------------------------------------------------
486! Surface downwelling longwave radiation.
487!-----------------------------------------------------------------------
488!
489 CALL get_2dfld (ng, itlm, idldwn, frcncid(idldwn,ng), &
490# if defined PIO_LIB && defined DISTRIBUTE
491 & frcpiofile(idldwn,ng), &
492# endif
493 & nffiles(ng), frc(1,ng), update(1), &
494 & lbi, ubi, lbj, ubj, 2, 1, &
495# ifdef MASKING
496 & grid(ng) % rmask, &
497# endif
498 & forces(ng) % lrflxG)
499 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
500# endif
501
502# if !defined ANA_TAIR && \
503 ((defined bulk_fluxes && !defined FORWARD_FLUXES) || \
504 defined ecosim || \
505 (defined shortwave && defined ana_srflux && defined albedo))
506!
507!-----------------------------------------------------------------------
508! Surface air temperature.
509!-----------------------------------------------------------------------
510!
511 CALL get_2dfld (ng, itlm, idtair, frcncid(idtair,ng), &
512# if defined PIO_LIB && defined DISTRIBUTE
513 & frcpiofile(idtair,ng), &
514# endif
515 & nffiles(ng), frc(1,ng), update(1), &
516 & lbi, ubi, lbj, ubj, 2, 1, &
517# ifdef MASKING
518 & grid(ng) % rmask, &
519# endif
520 & forces(ng) % TairG)
521 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
522# endif
523
524# if !defined ANA_HUMIDITY && \
525 ((defined bulk_fluxes && !defined FORWARD_FLUXES) || \
526 defined ecosim)
527!
528!-----------------------------------------------------------------------
529! Surface air humidity.
530!-----------------------------------------------------------------------
531!
532 CALL get_2dfld (ng, itlm, idqair, frcncid(idqair,ng), &
533# if defined PIO_LIB && defined DISTRIBUTE
534 & frcpiofile(idqair,ng), &
535# endif
536 & nffiles(ng), frc(1,ng), update(1), &
537 & lbi, ubi, lbj, ubj, 2, 1, &
538# ifdef MASKING
539 & grid(ng) % rmask, &
540# endif
541 & forces(ng) % HairG)
542 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
543# endif
544
545# if !defined ANA_RAIN && \
546 (defined bulk_fluxes && !defined FORWARD_FLUXES)
547!
548!-----------------------------------------------------------------------
549! Rain fall rate.
550!-----------------------------------------------------------------------
551!
552 CALL get_2dfld (ng, itlm, idrain, frcncid(idrain,ng), &
553# if defined PIO_LIB && defined DISTRIBUTE
554 & frcpiofile(idrain,ng), &
555# endif
556 & nffiles(ng), frc(1,ng), update(1), &
557 & lbi, ubi, lbj, ubj, 2, 1, &
558# ifdef MASKING
559 & grid(ng) % rmask, &
560# endif
561 & forces(ng) % rainG)
562 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
563# endif
564
565# ifndef ANA_STFLUX
566# if !(defined BULK_FLUXES || defined FORWARD_FLUXES)
567!
568!-----------------------------------------------------------------------
569! Surface net heat flux from input FRC file.
570!-----------------------------------------------------------------------
571!
572 CALL get_2dfld (ng, itlm, idtsur(itemp), &
573 & frcncid(idtsur(itemp),ng), &
574# if defined PIO_LIB && defined DISTRIBUTE
575 & frcpiofile(idtsur(itemp),ng), &
576# endif
577 & nffiles(ng), frc(1,ng), update(1), &
578 & lbi, ubi, lbj, ubj, 2, 1, &
579# ifdef MASKING
580 & grid(ng) % rmask, &
581# endif
582 & forces(ng) % stfluxG(:,:,:,itemp))
583 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
584
585# elif defined FORWARD_FLUXES
586!
587!-----------------------------------------------------------------------
588! Surface net heat flux from NLM forward file.
589!-----------------------------------------------------------------------
590!
591 CALL get_2dfld (ng, itlm, idtsur(itemp), blk(ng)%ncid, &
592# if defined PIO_LIB && defined DISTRIBUTE
593 & blk(ng)%pioFile, &
594# endif
595 & 1, blk(ng), update(1), &
596 & lbi, ubi, lbj, ubj, 2, 1, &
597# ifdef MASKING
598 & grid(ng) % rmask, &
599# endif
600 & forces(ng) % stfluxG(:,:,:,itemp))
601 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
602# endif
603# endif
604
605# if !defined ANA_SST && defined QCORRECTION
606!
607!-----------------------------------------------------------------------
608! Surface net heat flux correction fields: sea surface temperature
609! (SST).
610!-----------------------------------------------------------------------
611!
612 CALL get_2dfld (ng, itlm, idsstc, frcncid(idsstc,ng), &
613# if defined PIO_LIB && defined DISTRIBUTE
614 & frcpiofile(idsstc,ng), &
615# endif
616 & nffiles(ng), frc(1,ng), update(1), &
617 & lbi, ubi, lbj, ubj, 2, 1, &
618# ifdef MASKING
619 & grid(ng) % rmask, &
620# endif
621 & forces(ng) % sstG)
622 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
623# endif
624
625# if !defined ANA_DQDSST && defined QCORRECTION
626!
627!-----------------------------------------------------------------------
628! Surface net heat flux correction fields: heat flux sensitivity to
629! SST (dQdSST).
630!-----------------------------------------------------------------------
631!
632
633 CALL get_2dfld (ng, itlm, iddqdt, frcncid(iddqdt,ng), &
634# if defined PIO_LIB && defined DISTRIBUTE
635 & frcpiofile(iddqdt,ng), &
636# endif
637 & nffiles(ng), frc(1,ng), update(1), &
638 & lbi, ubi, lbj, ubj, 2, 1, &
639# ifdef MASKING
640 & grid(ng) % rmask, &
641# endif
642 & forces(ng) % dqdtG)
643 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
644# endif
645
646# ifndef ANA_BTFLUX
647!
648!-----------------------------------------------------------------------
649! Bottom net heat flux.
650!-----------------------------------------------------------------------
651!
652 CALL get_2dfld (ng, itlm, idtbot(itemp), &
653 & frcncid(idtbot(itemp),ng), &
654# if defined PIO_LIB && defined DISTRIBUTE
655 & frcpiofile(idtbot(itemp),ng), &
656# endif
657 & nffiles(ng), frc(1,ng), update(1), &
658 & lbi, ubi, lbj, ubj, 2, 1, &
659# ifdef MASKING
660 & grid(ng) % rmask, &
661# endif
662 & forces(ng) % btfluxG(:,:,:,itemp))
663 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
664# endif
665
666# if defined SALINITY && !defined ANA_SSFLUX
667# if !(defined EMINUSP || defined FORWARD_FLUXES || \
668 defined frc_coupling || defined srelaxation)
669!
670!-----------------------------------------------------------------------
671! Surface net freshwater flux: E-P from FRC file.
672!-----------------------------------------------------------------------
673!
674 CALL get_2dfld (ng, itlm, idsfwf, frcncid(idsfwf,ng), &
675# if defined PIO_LIB && defined DISTRIBUTE
676 & frcpiofile(idsfwf,ng), &
677# endif
678 & nffiles(ng), frc(1,ng), update(1), &
679 & lbi, ubi, lbj, ubj, 2, 1, &
680# ifdef MASKING
681 & grid(ng) % rmask, &
682# endif
683 & forces(ng) % stfluxG(:,:,:,isalt))
684 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
685
686# elif (defined EMINUSP || defined FORWARD_FLUXES || \
687 defined frc_coupling)
688!
689!-----------------------------------------------------------------------
690! Surface net freshwater flux (E-P) from NLM forward file.
691!-----------------------------------------------------------------------
692!
693 CALL get_2dfld (ng, itlm, idempf, blk(ng)%ncid, &
694# if defined PIO_LIB && defined DISTRIBUTE
695 & blk(ng)%pioFile, &
696# endif
697 & 1, blk(ng), update(1), &
698 & lbi, ubi, lbj, ubj, 2, 1, &
699# ifdef MASKING
700 & grid(ng) % rmask, &
701# endif
702 & forces(ng) % stfluxG(:,:,:,isalt))
703 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
704# endif
705
706# if !defined ANA_SSS && (defined SCORRECTION || defined SRELAXATION)
707!
708!-----------------------------------------------------------------------
709! Surface net freshwater flux correction field: sea surface salinity.
710!-----------------------------------------------------------------------
711!
712 CALL get_2dfld (ng, itlm, idsssc, frcncid(idsssc,ng), &
713# if defined PIO_LIB && defined DISTRIBUTE
714 & frcpiofile(idsssc,ng), &
715# endif
716 & nffiles(ng), frc(1,ng), update(1), &
717 & lbi, ubi, lbj, ubj, 2, 1, &
718# ifdef MASKING
719 & grid(ng) % rmask, &
720# endif
721 & forces(ng) % sssG)
722 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
723# endif
724
725# ifndef ANA_BSFLUX
726!
727!-----------------------------------------------------------------------
728! Bottom net freshwater flux.
729!-----------------------------------------------------------------------
730!
731 CALL get_2dfld (ng, itlm, idtbot(isalt), &
732 & frcncid(idtbot(isalt),ng), &
733# if defined PIO_LIB && defined DISTRIBUTE
734 & frcpiofile(idtbot(isalt),ng), &
735# endif
736 & nffiles(ng), frc(1,ng), update(1), &
737 & lbi, ubi, lbj, ubj, 2, 1, &
738# ifdef MASKING
739 & grid(ng) % rmask, &
740# endif
741 & forces(ng) % btfluxG(:,:,:,isalt))
742 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
743# endif
744# endif
745
746# if defined BIOLOGY || defined SEDIMENT || defined T_PASSIVE
747# ifndef ANA_SPFLUX
748!
749!-----------------------------------------------------------------------
750! Passive tracers surface fluxes.
751!-----------------------------------------------------------------------
752!
753 DO i=nat+1,nt(ng)
754 CALL get_2dfld (ng, itlm, idtsur(i), frcncid(idtsur(i),ng), &
755# if defined PIO_LIB && defined DISTRIBUTE
756 & frcpiofile(idtsur(i),ng), &
757# endif
758 & nffiles(ng), frc(1,ng), update(1), &
759 & lbi, ubi, lbj, ubj, 2, 1, &
760# ifdef MASKING
761 & grid(ng) % rmask, &
762# endif
763 & forces(ng) % stfluxG(:,:,:,i))
764 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
765 END DO
766# endif
767
768# ifndef ANA_BPFLUX
769!
770!-----------------------------------------------------------------------
771! Passive tracers bottom fluxes.
772!-----------------------------------------------------------------------
773!
774 DO i=nat+1,nt(ng)
775 CALL get_2dfld (ng, itlm, idtbot(i), frcncid(idtbot(i),ng), &
776# if defined PIO_LIB && defined DISTRIBUTE
777 & frcpiofile(idtbot(i),ng), &
778# endif
779 & nffiles(ng), frc(1,ng), update(1), &
780 & lbi, ubi, lbj, ubj, 2, 1, &
781# ifdef MASKING
782 & grid(ng) % rmask, &
783# endif
784 & forces(ng) % btfluxG(:,:,:,i))
785 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
786 END DO
787# endif
788# endif
789# endif
790!
791!=======================================================================
792! Read in open boundary conditions from BOUNDARY NetCDF file.
793!=======================================================================
794
795# ifndef ANA_FSOBC
796!
797 IF (lprocessobc(ng)) THEN
798 IF (tl_lbc(iwest,isfsur,ng)%acquire) THEN
799 CALL get_ngfld (ng, itlm, idzbry(iwest), &
800 & bryncid(idzbry(iwest),ng), &
801# if defined PIO_LIB && defined DISTRIBUTE
802 & brypiofile(idzbry(iwest),ng), &
803# endif
804 & nbcfiles(ng), bry(1,ng), &
805 & recordless, update(1), &
806 & jlb, jub, 1, 2, 0, mm(ng)+1, 1, &
807 & boundary(ng) % zetaG_west)
808 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
809 END IF
810!
811 IF (tl_lbc(ieast,isfsur,ng)%acquire) THEN
812 CALL get_ngfld (ng, itlm, idzbry(ieast), &
813 & bryncid(idzbry(ieast),ng), &
814# if defined PIO_LIB && defined DISTRIBUTE
815 & brypiofile(idzbry(ieast),ng), &
816# endif
817 & nbcfiles(ng), bry(1,ng), &
818 & recordless, update(1), &
819 & jlb, jub, 1, 2, 0, mm(ng)+1, 1, &
820 & boundary(ng) % zetaG_east)
821 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
822 END IF
823!
824 IF (tl_lbc(isouth,isfsur,ng)%acquire) THEN
825 CALL get_ngfld (ng, itlm, idzbry(isouth), &
826 & bryncid(idzbry(isouth),ng), &
827# if defined PIO_LIB && defined DISTRIBUTE
828 & brypiofile(idzbry(isouth),ng), &
829# endif
830 & nbcfiles(ng), bry(1,ng), &
831 & recordless, update(1), &
832 & ilb, iub, 1, 2, 0, lm(ng)+1, 1, &
833 & boundary(ng) % zetaG_south)
834 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
835 END IF
836!
837 IF (tl_lbc(inorth,isfsur,ng)%acquire) THEN
838 CALL get_ngfld (ng, itlm, idzbry(inorth), &
839 & bryncid(idzbry(inorth),ng), &
840# if defined PIO_LIB && defined DISTRIBUTE
841 & brypiofile(idzbry(inorth),ng), &
842# endif
843 & nbcfiles(ng), bry(1,ng), &
844 & recordless, update(1), &
845 & ilb, iub, 1, 2, 0, lm(ng)+1, 1, &
846 & boundary(ng) % zetaG_north)
847 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
848 END IF
849 END IF
850# endif
851
852# ifndef ANA_M2OBC
853!
854 IF (lprocessobc(ng)) THEN
855 IF (tl_lbc(iwest,isubar,ng)%acquire) THEN
856 CALL get_ngfld (ng, itlm, idu2bc(iwest), &
857 & bryncid(idu2bc(iwest),ng), &
858# if defined PIO_LIB && defined DISTRIBUTE
859 & brypiofile(idu2bc(iwest),ng), &
860# endif
861 & nbcfiles(ng), bry(1,ng), &
862 & recordless, update(1), &
863 & jlb, jub, 1, 2, 0, mm(ng)+1, 1, &
864 & boundary(ng) % ubarG_west)
865 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
866 END IF
867!
868 IF (tl_lbc(iwest,isvbar,ng)%acquire) THEN
869 CALL get_ngfld (ng, itlm, idv2bc(iwest), &
870 & bryncid(idv2bc(iwest),ng), &
871# if defined PIO_LIB && defined DISTRIBUTE
872 & brypiofile(idv2bc(iwest),ng), &
873# endif
874 & nbcfiles(ng), bry(1,ng), &
875 & recordless, update(1), &
876 & jlb, jub, 1, 2, 1, mm(ng)+1, 1, &
877 & boundary(ng) % vbarG_west)
878 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
879 END IF
880!
881 IF (tl_lbc(ieast,isubar,ng)%acquire) THEN
882 CALL get_ngfld (ng, itlm, idu2bc(ieast), &
883 & bryncid(idu2bc(ieast),ng), &
884# if defined PIO_LIB && defined DISTRIBUTE
885 & brypiofile(idu2bc(ieast),ng), &
886# endif
887 & nbcfiles(ng), bry(1,ng), &
888 & recordless, update(1), &
889 & jlb, jub, 1, 2, 0, mm(ng)+1, 1, &
890 & boundary(ng) % ubarG_east)
891 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
892 END IF
893!
894 IF (tl_lbc(ieast,isvbar,ng)%acquire) THEN
895 CALL get_ngfld (ng, itlm, idv2bc(ieast), &
896 & bryncid(idv2bc(ieast),ng), &
897# if defined PIO_LIB && defined DISTRIBUTE
898 & brypiofile(idv2bc(ieast),ng), &
899# endif
900 & nbcfiles(ng), bry(1,ng), &
901 & recordless, update(1), &
902 & jlb, jub, 1, 2, 1, mm(ng)+1, 1, &
903 & boundary(ng) % vbarG_east)
904 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
905 END IF
906!
907 IF (tl_lbc(isouth,isubar,ng)%acquire) THEN
908 CALL get_ngfld (ng, itlm, idu2bc(isouth), &
909 & bryncid(idu2bc(isouth),ng), &
910# if defined PIO_LIB && defined DISTRIBUTE
911 & brypiofile(idu2bc(isouth),ng), &
912# endif
913 & nbcfiles(ng), bry(1,ng), &
914 & recordless, update(1), &
915 & ilb, iub, 1, 2, 1, lm(ng)+1, 1, &
916 & boundary(ng) % ubarG_south)
917 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
918 END IF
919!
920 IF (tl_lbc(isouth,isvbar,ng)%acquire) THEN
921 CALL get_ngfld (ng, itlm, idv2bc(isouth), &
922 & bryncid(idv2bc(isouth),ng), &
923# if defined PIO_LIB && defined DISTRIBUTE
924 & brypiofile(idv2bc(isouth),ng), &
925# endif
926 & nbcfiles(ng), bry(1,ng), &
927 & recordless, update(1), &
928 & ilb, iub, 1, 2, 0, lm(ng)+1, 1, &
929 & boundary(ng) % vbarG_south)
930 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
931 END IF
932!
933 IF (tl_lbc(inorth,isubar,ng)%acquire) THEN
934 CALL get_ngfld (ng, itlm, idu2bc(inorth), &
935 & bryncid(idu2bc(inorth),ng), &
936# if defined PIO_LIB && defined DISTRIBUTE
937 & brypiofile(idu2bc(inorth),ng), &
938# endif
939 & nbcfiles(ng), bry(1,ng), &
940 & recordless, update(1), &
941 & ilb, iub, 1, 2, 1, lm(ng)+1, 1, &
942 & boundary(ng) % ubarG_north)
943 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
944 END IF
945!
946 IF (tl_lbc(inorth,isvbar,ng)%acquire) THEN
947 CALL get_ngfld (ng, itlm, idv2bc(inorth), &
948 & bryncid(idv2bc(inorth),ng), &
949# if defined PIO_LIB && defined DISTRIBUTE
950 & brypiofile(idv2bc(inorth),ng), &
951# endif
952 & nbcfiles(ng), bry(1,ng), &
953 & recordless, update(1), &
954 & ilb, iub, 1, 2, 0, lm(ng)+1, 1, &
955 & boundary(ng) % vbarG_north)
956 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
957 END IF
958 END IF
959# endif
960
961# ifdef SOLVE3D
962# ifndef ANA_M3OBC
963!
964 IF (lprocessobc(ng)) THEN
965 IF (tl_lbc(iwest,isuvel,ng)%acquire) THEN
966 CALL get_ngfld (ng, itlm, idu3bc(iwest), &
967 & bryncid(idu3bc(iwest),ng), &
968# if defined PIO_LIB && defined DISTRIBUTE
969 & brypiofile(idu3bc(iwest),ng), &
970# endif
971 & nbcfiles(ng), bry(1,ng), &
972 & recordless, update(1), &
973 & jlb, jub, n(ng), 2, 0, mm(ng)+1, n(ng), &
974 & boundary(ng) % uG_west)
975 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
976 END IF
977!
978 IF (tl_lbc(iwest,isvvel,ng)%acquire) THEN
979 CALL get_ngfld (ng, itlm, idv3bc(iwest), &
980 & bryncid(idv3bc(iwest),ng), &
981# if defined PIO_LIB && defined DISTRIBUTE
982 & brypiofile(idv3bc(iwest),ng), &
983# endif
984 & nbcfiles(ng), bry(1,ng), &
985 & recordless, update(1), &
986 & jlb, jub, n(ng), 2, 1, mm(ng)+1, n(ng), &
987 & boundary(ng) % vG_west)
988 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
989 END IF
990!
991 IF (tl_lbc(ieast,isuvel,ng)%acquire) THEN
992 CALL get_ngfld (ng, itlm, idu3bc(ieast), &
993 & bryncid(idu3bc(ieast),ng), &
994# if defined PIO_LIB && defined DISTRIBUTE
995 & brypiofile(idu3bc(ieast),ng), &
996# endif
997 & nbcfiles(ng), bry(1,ng), &
998 & recordless, update(1), &
999 & jlb, jub, n(ng), 2, 0, mm(ng)+1, n(ng), &
1000 & boundary(ng) % uG_east)
1001 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1002 END IF
1003!
1004 IF (tl_lbc(ieast,isvvel,ng)%acquire) THEN
1005 CALL get_ngfld (ng, itlm, idv3bc(ieast), &
1006 & bryncid(idv3bc(ieast),ng), &
1007# if defined PIO_LIB && defined DISTRIBUTE
1008 & brypiofile(idv3bc(ieast),ng), &
1009# endif
1010 & nbcfiles(ng), bry(1,ng), &
1011 & recordless, update(1), &
1012 & jlb, jub, n(ng), 2, 1, mm(ng)+1, n(ng), &
1013 & boundary(ng) % vG_east)
1014 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1015 END IF
1016!
1017 IF (tl_lbc(isouth,isuvel,ng)%acquire) THEN
1018 CALL get_ngfld (ng, itlm, idu3bc(isouth), &
1019 & bryncid(idu3bc(isouth),ng), &
1020# if defined PIO_LIB && defined DISTRIBUTE
1021 & brypiofile(idu3bc(isouth),ng), &
1022# endif
1023 & nbcfiles(ng), bry(1,ng), &
1024 & recordless, update(1), &
1025 & ilb, iub, n(ng), 2, 1, lm(ng)+1, n(ng), &
1026 & boundary(ng) % uG_south)
1027 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1028 END IF
1029!
1030 IF (tl_lbc(isouth,isvvel,ng)%acquire) THEN
1031 CALL get_ngfld (ng, itlm, idv3bc(isouth), &
1032 & bryncid(idv3bc(isouth),ng), &
1033# if defined PIO_LIB && defined DISTRIBUTE
1034 & brypiofile(idv3bc(isouth),ng), &
1035# endif
1036 & nbcfiles(ng), bry(1,ng), &
1037 & recordless, update(1), &
1038 & ilb, iub, n(ng), 2, 0, lm(ng)+1, n(ng), &
1039 & boundary(ng) % vG_south)
1040 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1041 END IF
1042!
1043 IF (tl_lbc(inorth,isuvel,ng)%acquire) THEN
1044 CALL get_ngfld (ng, itlm, idu3bc(inorth), &
1045 & bryncid(idu3bc(inorth),ng), &
1046# if defined PIO_LIB && defined DISTRIBUTE
1047 & brypiofile(idu3bc(inorth),ng), &
1048# endif
1049 & nbcfiles(ng), bry(1,ng), &
1050 & recordless, update(1), &
1051 & ilb, iub, n(ng), 2, 1, lm(ng)+1, n(ng), &
1052 & boundary(ng) % uG_north)
1053 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1054 END IF
1055!
1056 IF (tl_lbc(inorth,isvvel,ng)%acquire) THEN
1057 CALL get_ngfld (ng, itlm, idv3bc(inorth), &
1058 & bryncid(idv3bc(inorth),ng), &
1059# if defined PIO_LIB && defined DISTRIBUTE
1060 & brypiofile(idv3bc(inorth),ng), &
1061# endif
1062 & nbcfiles(ng), bry(1,ng), &
1063 & recordless, update(1), &
1064 & ilb, iub, n(ng), 2, 0, lm(ng)+1, n(ng), &
1065 & boundary(ng) % vG_north)
1066 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1067 END IF
1068 END IF
1069# endif
1070
1071# ifndef ANA_TOBC
1072!
1073 IF (lprocessobc(ng)) THEN
1074 DO i=1,nt(ng)
1075 IF (tl_lbc(iwest,istvar(i),ng)%acquire) THEN
1076 CALL get_ngfld (ng, itlm, idtbry(iwest,i), &
1077 & bryncid(idtbry(iwest,i),ng), &
1078# if defined PIO_LIB && defined DISTRIBUTE
1079 & brypiofile(idtbry(iwest,i),ng), &
1080# endif
1081 & nbcfiles(ng), bry(1,ng), &
1082 & recordless, update(1), &
1083 & jlb, jub, n(ng), 2, 0, mm(ng)+1, n(ng), &
1084 & boundary(ng) % tG_west(:,:,:,i))
1085 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1086 END IF
1087 END DO
1088!
1089 DO i=1,nt(ng)
1090 IF (tl_lbc(ieast,istvar(i),ng)%acquire) THEN
1091 CALL get_ngfld (ng, itlm, idtbry(ieast,i), &
1092 & bryncid(idtbry(ieast,i),ng), &
1093# if defined PIO_LIB && defined DISTRIBUTE
1094 & brypiofile(idtbry(ieast,i),ng), &
1095# endif
1096 & nbcfiles(ng), bry(1,ng), &
1097 & recordless, update(1), &
1098 & jlb, jub, n(ng), 2, 0, mm(ng)+1, n(ng), &
1099 & boundary(ng) % tG_east(:,:,:,i))
1100 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1101 END IF
1102 END DO
1103!
1104 DO i=1,nt(ng)
1105 IF (tl_lbc(isouth,istvar(i),ng)%acquire) THEN
1106 CALL get_ngfld (ng, itlm, idtbry(isouth,i), &
1107 & bryncid(idtbry(isouth,i),ng), &
1108# if defined PIO_LIB && defined DISTRIBUTE
1109 & brypiofile(idtbry(isouth,i),ng), &
1110# endif
1111 & nbcfiles(ng), bry(1,ng), &
1112 & recordless, update(1), &
1113 & ilb, iub, n(ng), 2, 0, lm(ng)+1, n(ng), &
1114 & boundary(ng) % tG_south(:,:,:,i))
1115 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1116 END IF
1117 END DO
1118!
1119 DO i=1,nt(ng)
1120 IF (tl_lbc(inorth,istvar(i),ng)%acquire) THEN
1121 CALL get_ngfld (ng, itlm, idtbry(inorth,i), &
1122 & bryncid(idtbry(inorth,i),ng), &
1123# if defined PIO_LIB && defined DISTRIBUTE
1124 & brypiofile(idtbry(inorth,i),ng), &
1125# endif
1126 & nbcfiles(ng), bry(1,ng), &
1127 & recordless, update(1), &
1128 & ilb, iub, n(ng), 2, 0, lm(ng)+1, n(ng), &
1129 & boundary(ng) % tG_north(:,:,:,i))
1130 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1131 END IF
1132 END DO
1133 END IF
1134# endif
1135# endif
1136!
1137!=======================================================================
1138! Read in data from Climatology NetCDF file.
1139!=======================================================================
1140
1141# ifndef ANA_SSH
1142!
1143! Free-surface.
1144!
1145 IF (lsshclm(ng)) THEN
1146 CALL get_2dfld (ng, itlm, idsshc, clmncid(idsshc,ng), &
1147# if defined PIO_LIB && defined DISTRIBUTE
1148 & clmpiofile(idsshc,ng), &
1149# endif
1150 & nclmfiles(ng), clm(1,ng), update(1), &
1151 & lbi, ubi, lbj, ubj, 2, 1, &
1152# ifdef MASKING
1153 & grid(ng) % rmask, &
1154# endif
1155 & clima(ng) % sshG)
1156 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1157 END IF
1158# endif
1159# ifndef ANA_M2CLIMA
1160!
1161! 2D momentum.
1162!
1163 IF (lm2clm(ng)) THEN
1164 CALL get_2dfld (ng, itlm, idubcl, clmncid(idubcl,ng), &
1165# if defined PIO_LIB && defined DISTRIBUTE
1166 & clmpiofile(idubcl,ng), &
1167# endif
1168 & nclmfiles(ng), clm(1,ng), update(1), &
1169 & lbi, ubi, lbj, ubj, 2, 1, &
1170# ifdef MASKING
1171 & grid(ng) % umask, &
1172# endif
1173 & clima(ng) % ubarclmG)
1174 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1175!
1176 CALL get_2dfld (ng, itlm, idvbcl, clmncid(idvbcl,ng), &
1177# if defined PIO_LIB && defined DISTRIBUTE
1178 & clmpiofile(idvbcl,ng), &
1179# endif
1180 & nclmfiles(ng), clm(1,ng), update(1), &
1181 & lbi, ubi, lbj, ubj, 2, 1, &
1182# ifdef MASKING
1183 & grid(ng) % vmask, &
1184# endif
1185 & clima(ng) % vbarclmG)
1186 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1187 END IF
1188# endif
1189# ifdef SOLVE3D
1190# ifndef ANA_M3CLIMA
1191!
1192! 3D momentum.
1193!
1194 IF (lm3clm(ng)) THEN
1195 CALL get_3dfld (ng, itlm, iduclm, clmncid(iduclm,ng), &
1196# if defined PIO_LIB && defined DISTRIBUTE
1197 & clmpiofile(iduclm,ng), &
1198# endif
1199 & nclmfiles(ng), clm(1,ng), update(1), &
1200 & lbi, ubi, lbj, ubj, 1, n(ng), 2, 1, &
1201# ifdef MASKING
1202 & grid(ng) % umask, &
1203# endif
1204 & clima(ng) % uclmG)
1205 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1206!
1207 CALL get_3dfld (ng, itlm, idvclm, clmncid(idvclm,ng), &
1208# if defined PIO_LIB && defined DISTRIBUTE
1209 & clmpiofile(idvclm,ng), &
1210# endif
1211 & nclmfiles(ng), clm(1,ng), update(1), &
1212 & lbi, ubi, lbj, ubj, 1, n(ng), 2, 1, &
1213# ifdef MASKING
1214 & grid(ng) % vmask, &
1215# endif
1216 & clima(ng) % vclmG)
1217 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1218 END IF
1219# endif
1220# ifndef ANA_TCLIMA
1221!
1222! Tracers.
1223!
1224 ic=0
1225 DO i=1,nt(ng)
1226 IF (ltracerclm(i,ng)) THEN
1227 ic=ic+1
1228 CALL get_3dfld (ng, itlm, idtclm(i), &
1229 & clmncid(idtclm(i),ng), &
1230# if defined PIO_LIB && defined DISTRIBUTE
1231 & clmpiofile(idtclm(i),ng), &
1232# endif
1233 & nclmfiles(ng), clm(1,ng), update(1), &
1234 & lbi, ubi, lbj, ubj, 1, n(ng), 2, 1, &
1235# ifdef MASKING
1236 & grid(ng) % rmask, &
1237# endif
1238 & clima(ng) % tclmG(:,:,:,:,ic))
1239 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1240 END IF
1241 END DO
1242# endif
1243# endif
1244
1245# ifdef FORWARD_READ
1246!
1247!-----------------------------------------------------------------------
1248! Read in forward state solution.
1249!-----------------------------------------------------------------------
1250!
1251! Read in free-surface.
1252!
1253 CALL get_2dfld (ng, itlm, idfsur, fwd(ng)%ncid, &
1254# if defined PIO_LIB && defined DISTRIBUTE
1255 & fwd(ng)%pioFile, &
1256# endif
1257 & 1, fwd(ng), update(1), &
1258 & lbi, ubi, lbj, ubj, 2, 1, &
1259# ifdef MASKING
1260 & grid(ng) % rmask, &
1261# endif
1262 & ocean(ng) % zetaG)
1263 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1264!
1265! Read 2D momentum.
1266!
1267 CALL get_2dfld (ng, itlm, idubar, fwd(ng)%ncid, &
1268# if defined PIO_LIB && defined DISTRIBUTE
1269 & fwd(ng)%pioFile, &
1270# endif
1271 & 1, fwd(ng), update(1), &
1272 & lbi, ubi, lbj, ubj, 2, 1, &
1273# ifdef MASKING
1274 & grid(ng) % umask, &
1275# endif
1276 & ocean(ng) % ubarG)
1277 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1278
1279 CALL get_2dfld (ng, itlm, idvbar, fwd(ng)%ncid, &
1280# if defined PIO_LIB && defined DISTRIBUTE
1281 & fwd(ng)%pioFile, &
1282# endif
1283 & 1, fwd(ng), update(1), &
1284 & lbi, ubi, lbj, ubj, 2, 1, &
1285# ifdef MASKING
1286 & grid(ng) % vmask, &
1287# endif
1288 & ocean(ng) % vbarG)
1289 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1290
1291# ifdef FORWARD_RHS
1292!
1293! Read in variables associated with 2D right-hand-side terms.
1294!
1295 CALL get_2dfld (ng, itlm, idrzet, fwd(ng)%ncid, &
1296# if defined PIO_LIB && defined DISTRIBUTE
1297 & fwd(ng)%pioFile, &
1298# endif
1299 & 1, fwd(ng), update(1), &
1300 & lbi, ubi, lbj, ubj, 2, 1, &
1301# ifdef MASKING
1302 & grid(ng) % rmask, &
1303# endif
1304 & ocean(ng) % rzetaG)
1305 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1306!
1307 CALL get_2dfld (ng, itlm, idru2d, fwd(ng)%ncid, &
1308# if defined PIO_LIB && defined DISTRIBUTE
1309 & fwd(ng)%pioFile, &
1310# endif
1311 & 1, fwd(ng), update(1), &
1312 & lbi, ubi, lbj, ubj, 2, 1, &
1313# ifdef MASKING
1314 & grid(ng) % umask, &
1315# endif
1316 & ocean(ng) % rubarG)
1317 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1318!
1319 CALL get_2dfld (ng, itlm, idrv2d, fwd(ng)%ncid, &
1320# if defined PIO_LIB && defined DISTRIBUTE
1321 & fwd(ng)%pioFile, &
1322# endif
1323 & 1, fwd(ng), update(1), &
1324 & lbi, ubi, lbj, ubj, 2, 1, &
1325# ifdef MASKING
1326 & grid(ng) % vmask, &
1327# endif
1328 & ocean(ng) % rvbarG)
1329 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1330# endif
1331
1332# ifdef SOLVE3D
1333!
1334! Read in variables associated with time-averaged 2D momentum terms.
1335!
1336 CALL get_2dfld (ng, itlm, idufx1, fwd(ng)%ncid, &
1337# if defined PIO_LIB && defined DISTRIBUTE
1338 & fwd(ng)%pioFile, &
1339# endif
1340 & 1, fwd(ng), update(1), &
1341 & lbi, ubi, lbj, ubj, 2, 1, &
1342# ifdef MASKING
1343 & grid(ng) % umask, &
1344# endif
1345 & coupling(ng) % DU_avg1G)
1346 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1347!
1348 CALL get_2dfld (ng, itlm, idufx2, fwd(ng)%ncid, &
1349# if defined PIO_LIB && defined DISTRIBUTE
1350 & fwd(ng)%pioFile, &
1351# endif
1352 & 1, fwd(ng), update(1), &
1353 & lbi, ubi, lbj, ubj, 2, 1, &
1354# ifdef MASKING
1355 & grid(ng) % umask, &
1356# endif
1357 & coupling(ng) % DU_avg2G)
1358 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1359!
1360 CALL get_2dfld (ng, itlm, idvfx1, fwd(ng)%ncid, &
1361# if defined PIO_LIB && defined DISTRIBUTE
1362 & fwd(ng)%pioFile, &
1363# endif
1364 & 1, fwd(ng), update(1), &
1365 & lbi, ubi, lbj, ubj, 2, 1, &
1366# ifdef MASKING
1367 & grid(ng) % vmask, &
1368# endif
1369 & coupling(ng) % DV_avg1G)
1370 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1371!
1372 CALL get_2dfld (ng, itlm, idvfx2, fwd(ng)%ncid, &
1373# if defined PIO_LIB && defined DISTRIBUTE
1374 & fwd(ng)%pioFile, &
1375# endif
1376 & 1, fwd(ng), update(1), &
1377 & lbi, ubi, lbj, ubj, 2, 1, &
1378# ifdef MASKING
1379 & grid(ng) % vmask, &
1380# endif
1381 & coupling(ng) % DV_avg2G)
1382 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1383!
1384! Read in 3D momentum.
1385!
1386 CALL get_3dfld (ng, itlm, iduvel, fwd(ng)%ncid, &
1387# if defined PIO_LIB && defined DISTRIBUTE
1388 & fwd(ng)%pioFile, &
1389# endif
1390 & 1, fwd(ng), update(1), &
1391 & lbi, ubi, lbj, ubj, 1, n(ng), 2, 1, &
1392# ifdef MASKING
1393 & grid(ng) % umask, &
1394# endif
1395 & ocean(ng) % uG)
1396 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1397!
1398 CALL get_3dfld (ng, irpm, idvvel, fwd(ng)%ncid, &
1399# if defined PIO_LIB && defined DISTRIBUTE
1400 & fwd(ng)%pioFile, &
1401# endif
1402 & 1, fwd(ng), update(1), &
1403 & lbi, ubi, lbj, ubj, 1, n(ng), 2, 1, &
1404# ifdef MASKING
1405 & grid(ng) % vmask, &
1406# endif
1407 & ocean(ng) % vG)
1408 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1409
1410# ifdef FORWARD_RHS
1411!
1412! Read in variables associated with 3D momentum right-hand-side terms.
1413!
1414 CALL get_2dfld (ng, itlm, idruct, fwd(ng)%ncid, &
1415# if defined PIO_LIB && defined DISTRIBUTE
1416 & fwd(ng)%pioFile, &
1417# endif
1418 & 1, fwd(ng), update(1), &
1419 & lbi, ubi, lbj, ubj, 2, 1, &
1420# ifdef MASKING
1421 & grid(ng) % umask, &
1422# endif
1423 & coupling(ng) % rufrcG)
1424 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1425!
1426 CALL get_2dfld (ng, itlm, idrvct, fwd(ng)%ncid, &
1427# if defined PIO_LIB && defined DISTRIBUTE
1428 & fwd(ng)%pioFile, &
1429# endif
1430 & 1, fwd(ng), update(1), &
1431 & lbi, ubi, lbj, ubj, 2, 1, &
1432# ifdef MASKING
1433 & grid(ng) % vmask, &
1434# endif
1435 & coupling(ng) % rvfrcG)
1436 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1437!
1438 CALL get_3dfld (ng, itlm, idru3d, fwd(ng)%ncid, &
1439# if defined PIO_LIB && defined DISTRIBUTE
1440 & fwd(ng)%pioFile, &
1441# endif
1442 & 1, fwd(ng), update(1), &
1443 & lbi, ubi, lbj, ubj, 1, n(ng), 2, 1, &
1444# ifdef MASKING
1445 & grid(ng) % umask, &
1446# endif
1447 & ocean(ng) % ruG)
1448 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1449!
1450 CALL get_3dfld (ng, itlm, idrv3d, fwd(ng)%ncid, &
1451# if defined PIO_LIB && defined DISTRIBUTE
1452 & fwd(ng)%pioFile, &
1453# endif
1454 & 1, fwd(ng), update(1), &
1455 & lbi, ubi, lbj, ubj, 1, n(ng), 2, 1, &
1456# ifdef MASKING
1457 & grid(ng) % vmask, &
1458# endif
1459 & ocean(ng) % rvG)
1460 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1461# endif
1462!
1463! Read in 3D tracers.
1464!
1465 DO i=1,nt(ng)
1466 CALL get_3dfld (ng, itlm, idtvar(i), fwd(ng)%ncid, &
1467# if defined PIO_LIB && defined DISTRIBUTE
1468 & fwd(ng)%pioFile, &
1469# endif
1470 & 1, fwd(ng), update(1), &
1471 & lbi, ubi, lbj, ubj, 1, n(ng), 2, 1, &
1472# ifdef MASKING
1473 & grid(ng) % rmask, &
1474# endif
1475 & ocean(ng) % tG(:,:,:,:,i))
1476 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1477 END DO
1478
1479# ifdef FORWARD_MIXING
1480!
1481! Read in vertical mixing variables.
1482!
1483 DO i=1,nat
1484 scale=fscale(iddiff(i),ng) ! save and rescale
1485 fscale(iddiff(i),ng)=tl_akt_fac(i,ng)
1486 CALL get_3dfld (ng, itlm, iddiff(i), fwd(ng)%ncid, &
1487# if defined PIO_LIB && defined DISTRIBUTE
1488 & fwd(ng)%pioFile, &
1489# endif
1490 & 1, fwd(ng), update(1), &
1491 & lbi, ubi, lbj, ubj, 0, n(ng), 2, 1, &
1492# ifdef MASKING
1493 & grid(ng) % rmask, &
1494# endif
1495 & mixing(ng) % AktG(:,:,:,:,i))
1496 fscale(iddiff(i),ng)=scale
1497 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1498 END DO
1499!
1500 scale=fscale(idvvis,ng) ! save and rescale
1501 fscale(idvvis,ng)=tl_akv_fac(ng)
1502 CALL get_3dfld (ng, itlm, idvvis, fwd(ng)%ncid, &
1503# if defined PIO_LIB && defined DISTRIBUTE
1504 & fwd(ng)%pioFile, &
1505# endif
1506 & 1, fwd(ng), update(1), &
1507 & lbi, ubi, lbj, ubj, 0, n(ng), 2, 1, &
1508# ifdef MASKING
1509 & grid(ng) % rmask, &
1510# endif
1511 & mixing(ng) % AkvG)
1512 fscale(idvvis,ng)=scale
1513 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1514# endif
1515
1516# if defined MY25_MIXING_NOT_YET || defined GLS_MIXING_NOT_YET
1517!
1518! Read in turbulent kinetic energy.
1519!
1520 CALL get_3dfld (ng, itlm, idmtke, fwd(ng)%ncid, &
1521# if defined PIO_LIB && defined DISTRIBUTE
1522 & fwd(ng)%pioFile, &
1523# endif
1524 & 1, fwd(ng), update(1), &
1525 & lbi, ubi, lbj, ubj, 0, n(ng), 2, 1, &
1526# ifdef MASKING
1527 & grid(ng) % rmask, &
1528# endif
1529 & mixing(ng) % tkeG)
1530 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1531!
1532! Read in turbulent kinetic energy times length scale.
1533!
1534 CALL get_3dfld (ng, itlm, idmtls, fwd(ng)%ncid, &
1535# if defined PIO_LIB && defined DISTRIBUTE
1536 & fwd(ng)%pioFile, &
1537# endif
1538 & 1, fwd(ng), update(1), &
1539 & lbi, ubi, lbj, ubj, 0, n(ng), 2, 1, &
1540# ifdef MASKING
1541 & grid(ng) % rmask, &
1542# endif
1543 & mixing(ng) % glsG)
1544 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1545!
1546! Read in vertical mixing length scale.
1547!
1548 CALL get_3dfld (ng, itlm, idvmls, fwd(ng)%ncid, &
1549# if defined PIO_LIB && defined DISTRIBUTE
1550 & fwd(ng)%pioFile, &
1551# endif
1552 & 1, fwd(ng), update(1), &
1553 & lbi, ubi, lbj, ubj, 0, n(ng), 2, 1, &
1554# ifdef MASKING
1555 & grid(ng) % rmask, &
1556# endif
1557 & mixing(ng) % LscaleG)
1558 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1559!
1560! Read in vertical mixing coefficient for turbulent kinetic energy.
1561!
1562 CALL get_3dfld (ng, itlm, idvmkk, fwd(ng)%ncid, &
1563# if defined PIO_LIB && defined DISTRIBUTE
1564 & fwd(ng)%pioFile, &
1565# endif
1566 & 1, fwd(ng), update(1), &
1567 & lbi, ubi, lbj, ubj, 0, n(ng), 2, 1, &
1568# ifdef MASKING
1569 & grid(ng) % rmask, &
1570# endif
1571 & mixing(ng) % AkkG)
1572 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1573
1574# ifdef GLS_MIXING_NOT_YET
1575!
1576! Read in vertical mixing coefficient for turbulent length scale.
1577!
1578 CALL get_3dfld (ng, itlm, idvmkp, fwd(ng)%ncid, &
1579# if defined PIO_LIB && defined DISTRIBUTE
1580 & fwd(ng)%pioFile, &
1581# endif
1582 & 1, fwd(ng), update(1), &
1583 & lbi, ubi, lbj, ubj, 0, n(ng), 2, 1, &
1584# ifdef MASKING
1585 & grid(ng) % rmask, &
1586# endif
1587 & mixing(ng) % AkpG)
1588 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1589# endif
1590# endif
1591
1592# ifdef LMD_MIXING_NOT_YET
1593!
1594! Read in depth of surface oceanic boundary layer.
1595!
1596 CALL get_2dfld (ng, itlm, idhsbl, fwd(ng)%ncid, &
1597# if defined PIO_LIB && defined DISTRIBUTE
1598 & fwd(ng)%pioFile, &
1599# endif
1600 & 1, fwd(ng), update(1), &
1601 & lbi, ubi, lbj, ubj, 2, 1, &
1602# ifdef MASKING
1603 & grid(ng) % rmask, &
1604# endif
1605 & mixing(ng) % hsblG)
1606 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1607# endif
1608
1609# ifdef LMD_BKPP_NOT_YET
1610!
1611! Read in depth of bottom oceanic boundary layer.
1612!
1613 CALL get_2dfld (ng, itlm, idhbbl, fwd(ng)%ncid, &
1614# if defined PIO_LIB && defined DISTRIBUTE
1615 & fwd(ng)%pioFile, &
1616# endif
1617 & 1, fwd(ng), update(1), &
1618 & lbi, ubi, lbj, ubj, 2, 1, &
1619# ifdef MASKING
1620 & grid(ng) % rmask, &
1621# endif
1622 & mixing(ng) % hbblG)
1623 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1624# endif
1625
1626# ifdef LMD_NONLOCAL_NOT_YET
1627!
1628! Read in boundary layer nonlocal transport.
1629!
1630 DO i=1,nat
1631 CALL get_3dfld (ng, itlm, idghat(i), fwd(ng)%ncid, &
1632# if defined PIO_LIB && defined DISTRIBUTE
1633 & fwd(ng)%pioFile, &
1634# endif
1635 & 1, fwd(ng), update(1), &
1636 & lbi, ubi, lbj, ubj, 0, n(ng), 2, 1, &
1637# ifdef MASKING
1638 & grid(ng) % rmask, &
1639# endif
1640 & mixing(ng) % ghatsG(:,:,:,i))
1641 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1642 END DO
1643# endif
1644# endif
1645
1646# ifdef WEAK_CONSTRAINT
1647!
1648!-----------------------------------------------------------------------
1649! Read frequent impulse forcing for weak constraint.
1650!-----------------------------------------------------------------------
1651!
1652 IF (frequentimpulse(ng)) THEN
1653 CALL get_2dfld (ng, itlm, idztlf, tlf(ng)%ncid, &
1654# if defined PIO_LIB && defined DISTRIBUTE
1655 & tlf(ng)%pioFile, &
1656# endif
1657 & 1, tlf(ng), update(1), &
1658 & lbi, ubi, lbj, ubj, 2, 1, &
1659# ifdef MASKING
1660 & grid(ng) % rmask, &
1661# endif
1662 & ocean(ng) % f_zetaG)
1663 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1664
1665# ifndef SOLVE3D
1666!
1667! Read 2D momentum forcing.
1668!
1669 CALL get_2dfld (ng, itlm, idubtf, tlf(ng)%ncid, &
1670# if defined PIO_LIB && defined DISTRIBUTE
1671 & tlf(ng)%pioFile, &
1672# endif
1673 & 1, tlf(ng), update(1), &
1674 & lbi, ubi, lbj, ubj, 2, 1, &
1675# ifdef MASKING
1676 & grid(ng) % umask, &
1677# endif
1678 & ocean(ng) % f_ubarG)
1679 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1680!
1681 CALL get_2dfld (ng, itlm, idvbtf, tlf(ng)%ncid, &
1682# if defined PIO_LIB && defined DISTRIBUTE
1683 & tlf(ng)%pioFile, &
1684# endif
1685 & 1, tlf(ng), update(1), &
1686 & lbi, ubi, lbj, ubj, 2, 1, &
1687# ifdef MASKING
1688 & grid(ng) % vmask, &
1689# endif
1690 & ocean(ng) % f_vbarG)
1691 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1692# endif
1693
1694# ifdef SOLVE3D
1695!
1696! Read in 3D momentum forcing.
1697!
1698 CALL get_3dfld (ng, itlm, idutlf, tlf(ng)%ncid, &
1699# if defined PIO_LIB && defined DISTRIBUTE
1700 & tlf(ng)%pioFile, &
1701# endif
1702 & 1, tlf(ng), update(1), &
1703 & lbi, ubi, lbj, ubj, 1, n(ng), 2, 1, &
1704# ifdef MASKING
1705 & grid(ng) % umask, &
1706# endif
1707 & ocean(ng) % f_uG)
1708 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1709!
1710 CALL get_3dfld (ng, itlm, idvtlf, tlf(ng)%ncid, &
1711# if defined PIO_LIB && defined DISTRIBUTE
1712 & tlf(ng)%pioFile, &
1713# endif
1714 & 1, tlf(ng), update(1), &
1715 & lbi, ubi, lbj, ubj, 1, n(ng), 2, 1, &
1716# ifdef MASKING
1717 & grid(ng) % vmask, &
1718# endif
1719 & ocean(ng) % f_vG)
1720 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1721!
1722! Read in 3D tracers forcing.
1723!
1724 DO i=1,nt(ng)
1725 CALL get_3dfld (ng, itlm, idttlf(i), tlf(ng)%ncid, &
1726# if defined PIO_LIB && defined DISTRIBUTE
1727 & tlf(ng)%pioFile, &
1728# endif
1729 & 1, tlf(ng), update(1), &
1730 & lbi, ubi, lbj, ubj, 1, n(ng), 2, 1, &
1731# ifdef MASKING
1732 & grid(ng) % rmask, &
1733# endif
1734 & ocean(ng) % f_tG(:,:,:,:,i))
1735 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1736 END DO
1737# endif
1738 END IF
1739# endif
1740# endif
1741
1742# ifdef PROFILE
1743!
1744!-----------------------------------------------------------------------
1745! Turn off input data time wall clock.
1746!-----------------------------------------------------------------------
1747!
1748 CALL wclock_off (ng, itlm, 3, __line__, myfile)
1749# endif
1750!
1751 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_3dfld(ng, model, ifield, ncid, piofile, nfiles, s, update, lbi, ubi, lbj, ubj, lbk, ubk, iout, irec, fmask, fout)
Definition get_3dfld.F:14
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
type(t_boundary), dimension(:), allocatable boundary
type(t_clima), dimension(:), allocatable clima
Definition mod_clima.F:153
type(t_coupling), dimension(:), allocatable coupling
type(t_forces), dimension(:), allocatable forces
Definition mod_forces.F:554
type(t_grid), dimension(:), allocatable grid
Definition mod_grid.F:365
type(t_io), dimension(:), allocatable ssf
type(t_io), dimension(:), allocatable tlf
integer, dimension(:), allocatable nclmfiles
type(t_io), dimension(:), allocatable blk
type(t_io), dimension(:,:), allocatable frc
type(t_io), dimension(:,:), allocatable bry
type(t_io), dimension(:), allocatable fwd
integer, dimension(:), allocatable nbcfiles
type(t_io), dimension(:,:), allocatable clm
integer, dimension(:), allocatable nffiles
type(t_mixing), dimension(:), allocatable mixing
Definition mod_mixing.F:399
integer iddqdt
integer idvair
integer idvmls
integer idubtf
integer, dimension(2) iddiff
integer, dimension(:,:), allocatable clmncid
type(file_desc_t), dimension(:,:), pointer frcpiofile
integer, dimension(4) idu3bc
integer idrtra
integer idwdis
integer idrv3d
integer idubar
integer, dimension(4) idzbry
integer idvvel
integer idhsbl
integer, dimension(4) idu2bc
integer idvsms
integer idcfra
integer idwlen
integer isvvel
integer isvbar
integer idvclm
integer, dimension(:), allocatable idrtrc
integer idpair
integer, dimension(:), allocatable idtbot
integer idvfx2
integer, dimension(:), allocatable idtsur
integer idru2d
integer idvmkp
integer idempf
integer, dimension(:), allocatable idtclm
integer, dimension(:,:), allocatable frcncid
integer idfsur
integer idwbrk
integer, dimension(:), allocatable idtvar
integer, dimension(:), allocatable istvar
integer idhbbl
integer, dimension(:,:), allocatable bryncid
integer idvfx1
integer idldwn
integer isuvel
integer idufx2
integer idsssc
integer isfsur
integer iduclm
integer idztlf
integer idsfwf
integer idvbtf
integer, dimension(4) idv3bc
real(dp), dimension(:,:), allocatable fscale
integer iduair
integer, dimension(:), allocatable idttlf
integer idmtke
integer iduvel
integer, dimension(2) idghat
integer idqair
integer idutlf
integer isubar
integer idubcl
integer, dimension(4) idv2bc
integer idlrad
integer idvtlf
integer idru3d
integer idusms
integer idvmkk
integer idwamp
integer idwdir
integer idvvis
integer idwptp
integer idrzet
type(file_desc_t), dimension(:,:), pointer clmpiofile
integer idrvct
integer idufx1
type(file_desc_t), dimension(:,:), pointer brypiofile
integer idrain
integer idvbcl
integer idsrad
integer idmtls
integer idsshc
integer idwpbt
integer idruct
integer idworb
integer idsstc
integer, dimension(:,:), allocatable idtbry
integer idtair
integer idrv2d
integer idvbar
type(t_ocean), dimension(:), allocatable ocean
Definition mod_ocean.F:351
integer nat
Definition mod_param.F:499
integer, parameter irpm
Definition mod_param.F:664
integer, dimension(:), allocatable n
Definition mod_param.F:479
type(t_bounds), dimension(:), allocatable bounds
Definition mod_param.F:232
type(t_lbc), dimension(:,:,:), allocatable tl_lbc
Definition mod_param.F:379
integer, dimension(:), allocatable lm
Definition mod_param.F:455
integer, parameter itlm
Definition mod_param.F:663
integer, dimension(:), allocatable nt
Definition mod_param.F:489
integer, dimension(:), allocatable mm
Definition mod_param.F:456
logical, dimension(:), allocatable luvsrc
logical, dimension(:,:), allocatable ltracersrc
logical, dimension(:), allocatable lprocessobc
integer, parameter iwest
logical, dimension(:), allocatable lm3clm
real(r8), dimension(:), allocatable tl_akv_fac
logical, dimension(:), allocatable lsshclm
logical, dimension(:), allocatable frequentimpulse
logical, dimension(:), allocatable lwsrc
integer exit_flag
integer isalt
integer itemp
integer, parameter isouth
logical, dimension(:), allocatable lm2clm
integer, parameter ieast
logical, dimension(:,:), allocatable ltracerclm
integer, parameter inorth
integer noerror
real(r8), dimension(:,:), allocatable tl_akt_fac
type(t_sources), dimension(:), allocatable sources
Definition mod_sources.F:90
integer, dimension(:), allocatable nsrc
Definition mod_sources.F:97
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_iounits::blk, mod_boundary::boundary, mod_param::bounds, mod_iounits::bry, mod_ncparam::bryncid, mod_ncparam::brypiofile, mod_clima::clima, mod_iounits::clm, mod_ncparam::clmncid, mod_ncparam::clmpiofile, mod_coupling::coupling, mod_scalars::exit_flag, mod_forces::forces, strings_mod::founderror(), mod_iounits::frc, mod_ncparam::frcncid, mod_ncparam::frcpiofile, mod_scalars::frequentimpulse, mod_ncparam::fscale, mod_iounits::fwd, get_2dfld(), get_3dfld(), get_ngfld(), mod_grid::grid, mod_ncparam::idcfra, mod_ncparam::iddiff, mod_ncparam::iddqdt, mod_ncparam::idempf, mod_ncparam::idfsur, mod_ncparam::idghat, mod_ncparam::idhbbl, mod_ncparam::idhsbl, mod_ncparam::idldwn, mod_ncparam::idlrad, mod_ncparam::idmtke, mod_ncparam::idmtls, mod_ncparam::idpair, mod_ncparam::idqair, mod_ncparam::idrain, mod_ncparam::idrtra, mod_ncparam::idrtrc, mod_ncparam::idru2d, mod_ncparam::idru3d, mod_ncparam::idruct, mod_ncparam::idrv2d, mod_ncparam::idrv3d, mod_ncparam::idrvct, mod_ncparam::idrzet, mod_ncparam::idsfwf, mod_ncparam::idsrad, mod_ncparam::idsshc, mod_ncparam::idsssc, mod_ncparam::idsstc, mod_ncparam::idtair, mod_ncparam::idtbot, mod_ncparam::idtbry, mod_ncparam::idtclm, mod_ncparam::idtsur, mod_ncparam::idttlf, mod_ncparam::idtvar, mod_ncparam::idu2bc, mod_ncparam::idu3bc, mod_ncparam::iduair, mod_ncparam::idubar, mod_ncparam::idubcl, mod_ncparam::idubtf, mod_ncparam::iduclm, mod_ncparam::idufx1, mod_ncparam::idufx2, mod_ncparam::idusms, mod_ncparam::idutlf, mod_ncparam::iduvel, mod_ncparam::idv2bc, mod_ncparam::idv3bc, mod_ncparam::idvair, mod_ncparam::idvbar, mod_ncparam::idvbcl, mod_ncparam::idvbtf, mod_ncparam::idvclm, mod_ncparam::idvfx1, mod_ncparam::idvfx2, mod_ncparam::idvmkk, mod_ncparam::idvmkp, mod_ncparam::idvmls, mod_ncparam::idvsms, mod_ncparam::idvtlf, mod_ncparam::idvvel, mod_ncparam::idvvis, mod_ncparam::idwamp, mod_ncparam::idwbrk, mod_ncparam::idwdir, mod_ncparam::idwdis, mod_ncparam::idwlen, mod_ncparam::idworb, mod_ncparam::idwpbt, mod_ncparam::idwptp, mod_ncparam::idzbry, mod_ncparam::idztlf, mod_scalars::ieast, mod_scalars::inorth, mod_param::irpm, mod_scalars::isalt, mod_ncparam::isfsur, mod_scalars::isouth, mod_ncparam::istvar, mod_ncparam::isubar, mod_ncparam::isuvel, mod_ncparam::isvbar, mod_ncparam::isvvel, mod_scalars::itemp, mod_param::itlm, mod_scalars::iwest, mod_param::lm, mod_scalars::lm2clm, mod_scalars::lm3clm, mod_scalars::lprocessobc, mod_scalars::lsshclm, mod_scalars::ltracerclm, mod_scalars::ltracersrc, mod_scalars::luvsrc, mod_scalars::lwsrc, mod_mixing::mixing, mod_param::mm, mod_param::n, mod_param::nat, mod_iounits::nbcfiles, mod_iounits::nclmfiles, mod_iounits::nffiles, mod_scalars::noerror, mod_sources::nsrc, mod_param::nt, mod_ocean::ocean, mod_sources::sources, mod_iounits::ssf, mod_scalars::tl_akt_fac, mod_scalars::tl_akv_fac, mod_param::tl_lbc, mod_iounits::tlf, wclock_off(), and wclock_on().

Referenced by propagator_mod::propagator_fsv(), propagator_mod::propagator_fte(), propagator_mod::propagator_hop(), propagator_mod::propagator_hso(), propagator_mod::propagator_op(), propagator_mod::propagator_so(), tl_initial(), tl_main3d(), and roms_kernel_mod::tlm_initial().

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