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

Go to the source code of this file.

Functions/Subroutines

subroutine rp_get_data (ng)
 

Function/Subroutine Documentation

◆ rp_get_data()

subroutine rp_get_data ( integer, intent(in) ng)

Definition at line 3 of file rp_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, irpm, 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, irpm, 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, irpm, 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, irpm, 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 , irpm, 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, irpm, 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, irpm, 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, irpm, 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, irpm, 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, irpm, 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 used in 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, irpm, 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, irpm, 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, irpm, 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, irpm, 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, irpm, 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, irpm, 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, irpm, 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, irpm, 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, irpm, 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, irpm, 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, irpm, 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, irpm, 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, irpm, 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, irpm, 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, irpm, 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, irpm, 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, irpm, 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, irpm, 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, irpm, 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 CALL get_2dfld (ng, irpm, iddqdt, frcncid(iddqdt,ng), &
633# if defined PIO_LIB && defined DISTRIBUTE
634 & frcpiofile(iddqdt,ng), &
635# endif
636 & nffiles(ng), frc(1,ng), update(1), &
637 & lbi, ubi, lbj, ubj, 2, 1, &
638# ifdef MASKING
639 & grid(ng) % rmask, &
640# endif
641 & forces(ng) % dqdtG)
642 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
643# endif
644
645# ifndef ANA_BTFLUX
646!
647!-----------------------------------------------------------------------
648! Bottom net heat flux.
649!-----------------------------------------------------------------------
650!
651 CALL get_2dfld (ng, irpm, idtbot(itemp), &
652 & frcncid(idtbot(itemp),ng), &
653# if defined PIO_LIB && defined DISTRIBUTE
654 & frcpiofile(idtbot(itemp),ng), &
655# endif
656 & nffiles(ng), frc(1,ng), update(1), &
657 & lbi, ubi, lbj, ubj, 2, 1, &
658# ifdef MASKING
659 & grid(ng) % rmask, &
660# endif
661 & forces(ng) % btfluxG(:,:,:,itemp))
662 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
663# endif
664
665# if defined SALINITY && !defined ANA_SSFLUX
666# if !(defined EMINUSP || defined FORWARD_FLUXES || \
667 defined frc_coupling || defined srelaxation)
668!
669!-----------------------------------------------------------------------
670! Surface net freshwater flux: E-P from FRC file.
671!-----------------------------------------------------------------------
672!
673 CALL get_2dfld (ng, irpm, idsfwf, frcncid(idsfwf,ng), &
674# if defined PIO_LIB && defined DISTRIBUTE
675 & frcpiofile(idsfwf,ng), &
676# endif
677 & nffiles(ng), frc(1,ng), update(1), &
678 & lbi, ubi, lbj, ubj, 2, 1, &
679# ifdef MASKING
680 & grid(ng) % rmask, &
681# endif
682 & forces(ng) % stfluxG(:,:,:,isalt))
683 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
684
685# elif (defined EMINUSP || defined FORWARD_FLUXES || \
686 defined frc_coupling)
687!
688!-----------------------------------------------------------------------
689! Surface net freshwater flux (E-P) from NLM forward file.
690!-----------------------------------------------------------------------
691!
692 CALL get_2dfld (ng, irpm, idempf, blk(ng)%ncid, &
693# if defined PIO_LIB && defined DISTRIBUTE
694 & blk(ng)%pioFile, &
695# endif
696 & 1, blk(ng), update(1), &
697 & lbi, ubi, lbj, ubj, 2, 1, &
698# ifdef MASKING
699 & grid(ng) % rmask, &
700# endif
701 & forces(ng) % stfluxG(:,:,:,isalt))
702 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
703# endif
704
705# if !defined ANA_SSS && (defined SCORRECTION || defined SRELAXATION)
706!
707!-----------------------------------------------------------------------
708! Surface net freshwater flux correction field: sea surface salinity.
709!-----------------------------------------------------------------------
710!
711 CALL get_2dfld (ng, irpm, idsssc, frcncid(idsssc,ng), &
712# if defined PIO_LIB && defined DISTRIBUTE
713 & frcpiofile(idsssc,ng), &
714# endif
715 & nffiles(ng), frc(1,ng), update(1), &
716 & lbi, ubi, lbj, ubj, 2, 1, &
717# ifdef MASKING
718 & grid(ng) % rmask, &
719# endif
720 & forces(ng) % sssG)
721 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
722# endif
723
724# ifndef ANA_BSFLUX
725!
726!-----------------------------------------------------------------------
727! Bottom net freshwater flux.
728!-----------------------------------------------------------------------
729!
730 CALL get_2dfld (ng, irpm, idtbot(isalt), &
731 & frcncid(idtbot(isalt),ng), &
732# if defined PIO_LIB && defined DISTRIBUTE
733 & frcpiofile(idtbot(isalt),ng), &
734# endif
735 & nffiles(ng), frc(1,ng), update(1), &
736 & lbi, ubi, lbj, ubj, 2, 1, &
737# ifdef MASKING
738 & grid(ng) % rmask, &
739# endif
740 & forces(ng) % btfluxG(:,:,:,isalt))
741 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
742# endif
743# endif
744
745# if defined BIOLOGY || defined SEDIMENT || defined T_PASSIVE
746# ifndef ANA_SPFLUX
747!
748!-----------------------------------------------------------------------
749! Passive tracers surface fluxes.
750!-----------------------------------------------------------------------
751!
752 DO i=nat+1,nt(ng)
753 CALL get_2dfld (ng, irpm, idtsur(i), frcncid(idtsur(i),ng), &
754# if defined PIO_LIB && defined DISTRIBUTE
755 & frcpiofile(idtsur(i),ng), &
756# endif
757 & nffiles(ng), frc(1,ng), update(1), &
758 & lbi, ubi, lbj, ubj, 2, 1, &
759# ifdef MASKING
760 & grid(ng) % rmask, &
761# endif
762 & forces(ng) % stfluxG(:,:,:,i))
763 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
764 END DO
765# endif
766
767# ifndef ANA_BPFLUX
768!
769!-----------------------------------------------------------------------
770! Passive tracers bottom fluxes.
771!-----------------------------------------------------------------------
772!
773 DO i=nat+1,nt(ng)
774 CALL get_2dfld (ng, irpm, idtbot(i), frcncid(idtbot(i),ng), &
775# if defined PIO_LIB && defined DISTRIBUTE
776 & frcpiofile(idtbot(i),ng), &
777# endif
778 & nffiles(ng), frc(1,ng), update(1), &
779 & lbi, ubi, lbj, ubj, 2, 1, &
780# ifdef MASKING
781 & grid(ng) % rmask, &
782# endif
783 & forces(ng) % btfluxG(:,:,:,i))
784 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
785 END DO
786# endif
787# endif
788# endif
789!
790!=======================================================================
791! Read in open boundary conditions from BOUNDARY NetCDF file.
792!=======================================================================
793
794# ifndef ANA_FSOBC
795!
796 IF (lprocessobc(ng)) THEN
797 IF (tl_lbc(iwest,isfsur,ng)%acquire) THEN
798 CALL get_ngfld (ng, irpm, idzbry(iwest), &
799 & bryncid(idzbry(iwest),ng), &
800# if defined PIO_LIB && defined DISTRIBUTE
801 & brypiofile(idzbry(iwest),ng), &
802# endif
803 & nbcfiles(ng), bry(1,ng), &
804 & recordless, update(1), &
805 & jlb, jub, 1, 2, 0, mm(ng)+1, 1, &
806 & boundary(ng) % zetaG_west)
807 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
808 END IF
809!
810 IF (tl_lbc(ieast,isfsur,ng)%acquire) THEN
811 CALL get_ngfld (ng, irpm, idzbry(ieast), &
812 & bryncid(idzbry(ieast),ng), &
813# if defined PIO_LIB && defined DISTRIBUTE
814 & brypiofile(idzbry(ieast),ng), &
815# endif
816 & nbcfiles(ng), bry(1,ng), &
817 & recordless, update(1), &
818 & jlb, jub, 1, 2, 0, mm(ng)+1, 1, &
819 & boundary(ng) % zetaG_east)
820 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
821 END IF
822!
823 IF (tl_lbc(isouth,isfsur,ng)%acquire) THEN
824 CALL get_ngfld (ng, irpm, idzbry(isouth), &
825 & bryncid(idzbry(isouth),ng), &
826# if defined PIO_LIB && defined DISTRIBUTE
827 & brypiofile(idzbry(isouth),ng), &
828# endif
829 & nbcfiles(ng), bry(1,ng), &
830 & recordless, update(1), &
831 & ilb, iub, 1, 2, 0, lm(ng)+1, 1, &
832 & boundary(ng) % zetaG_south)
833 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
834 END IF
835!
836 IF (tl_lbc(inorth,isfsur,ng)%acquire) THEN
837 CALL get_ngfld (ng, irpm, idzbry(inorth), &
838 & bryncid(idzbry(inorth),ng), &
839# if defined PIO_LIB && defined DISTRIBUTE
840 & brypiofile(idzbry(inorth),ng), &
841# endif
842 & nbcfiles(ng), bry(1,ng), &
843 & recordless, update(1), &
844 & ilb, iub, 1, 2, 0, lm(ng)+1, 1, &
845 & boundary(ng) % zetaG_north)
846 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
847 END IF
848 END IF
849# endif
850
851# ifndef ANA_M2OBC
852!
853 IF (lprocessobc(ng)) THEN
854 IF (tl_lbc(iwest,isubar,ng)%acquire) THEN
855 CALL get_ngfld (ng, irpm, idu2bc(iwest), &
856 & bryncid(idu2bc(iwest),ng), &
857# if defined PIO_LIB && defined DISTRIBUTE
858 & brypiofile(idu2bc(iwest),ng), &
859# endif
860 & nbcfiles(ng), bry(1,ng), &
861 & recordless, update(1), &
862 & jlb, jub, 1, 2, 0, mm(ng)+1, 1, &
863 & boundary(ng) % ubarG_west)
864 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
865 END IF
866!
867 IF (tl_lbc(iwest,isvbar,ng)%acquire) THEN
868 CALL get_ngfld (ng, irpm, idv2bc(iwest), &
869 & bryncid(idv2bc(iwest),ng), &
870# if defined PIO_LIB && defined DISTRIBUTE
871 & brypiofile(idv2bc(iwest),ng), &
872# endif
873 & nbcfiles(ng), bry(1,ng), &
874 & recordless, update(1), &
875 & jlb, jub, 1, 2, 1, mm(ng)+1, 1, &
876 & boundary(ng) % vbarG_west)
877 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
878 END IF
879!
880 IF (tl_lbc(ieast,isubar,ng)%acquire) THEN
881 CALL get_ngfld (ng, irpm, idu2bc(ieast), &
882 & bryncid(idu2bc(ieast),ng), &
883# if defined PIO_LIB && defined DISTRIBUTE
884 & brypiofile(idu2bc(ieast),ng), &
885# endif
886 & nbcfiles(ng), bry(1,ng), &
887 & recordless, update(1), &
888 & jlb, jub, 1, 2, 0, mm(ng)+1, 1, &
889 & boundary(ng) % ubarG_east)
890 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
891 END IF
892!
893 IF (tl_lbc(ieast,isvbar,ng)%acquire) THEN
894 CALL get_ngfld (ng, irpm, idv2bc(ieast), &
895 & bryncid(idv2bc(ieast),ng), &
896# if defined PIO_LIB && defined DISTRIBUTE
897 & brypiofile(idv2bc(ieast),ng), &
898# endif
899 & nbcfiles(ng), bry(1,ng), &
900 & recordless, update(1), &
901 & jlb, jub, 1, 2, 1, mm(ng)+1, 1, &
902 & boundary(ng) % vbarG_east)
903 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
904 END IF
905!
906 IF (tl_lbc(isouth,isubar,ng)%acquire) THEN
907 CALL get_ngfld (ng, irpm, idu2bc(isouth), &
908 & bryncid(idu2bc(isouth),ng), &
909# if defined PIO_LIB && defined DISTRIBUTE
910 & brypiofile(idu2bc(isouth),ng), &
911# endif
912 & nbcfiles(ng), bry(1,ng), &
913 & recordless, update(1), &
914 & ilb, iub, 1, 2, 1, lm(ng)+1, 1, &
915 & boundary(ng) % ubarG_south)
916 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
917 END IF
918!
919 IF (tl_lbc(isouth,isvbar,ng)%acquire) THEN
920 CALL get_ngfld (ng, irpm, idv2bc(isouth), &
921 & bryncid(idv2bc(isouth),ng), &
922# if defined PIO_LIB && defined DISTRIBUTE
923 & brypiofile(idv2bc(isouth),ng), &
924# endif
925 & nbcfiles(ng), bry(1,ng), &
926 & recordless, update(1), &
927 & ilb, iub, 1, 2, 0, lm(ng)+1, 1, &
928 & boundary(ng) % vbarG_south)
929 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
930 END IF
931!
932 IF (tl_lbc(inorth,isubar,ng)%acquire) THEN
933 CALL get_ngfld (ng, irpm, idu2bc(inorth), &
934 & bryncid(idu2bc(inorth),ng), &
935# if defined PIO_LIB && defined DISTRIBUTE
936 & brypiofile(idu2bc(inorth),ng), &
937# endif
938 & nbcfiles(ng), bry(1,ng), &
939 & recordless, update(1), &
940 & ilb, iub, 1, 2, 1, lm(ng)+1, 1, &
941 & boundary(ng) % ubarG_north)
942 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
943 END IF
944
945 IF (tl_lbc(inorth,isvbar,ng)%acquire) THEN
946 CALL get_ngfld (ng, irpm, idv2bc(inorth), &
947 & bryncid(idv2bc(inorth),ng), &
948# if defined PIO_LIB && defined DISTRIBUTE
949 & brypiofile(idv2bc(inorth),ng), &
950# endif
951 & nbcfiles(ng), bry(1,ng), &
952 & recordless, update(1), &
953 & ilb, iub, 1, 2, 0, lm(ng)+1, 1, &
954 & boundary(ng) % vbarG_north)
955 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
956 END IF
957 END IF
958# endif
959
960# ifdef SOLVE3D
961# ifndef ANA_M3OBC
962!
963 IF (lprocessobc(ng)) THEN
964 IF (tl_lbc(iwest,isuvel,ng)%acquire) THEN
965 CALL get_ngfld (ng, irpm, idu3bc(iwest), &
966 & bryncid(idu3bc(iwest),ng), &
967# if defined PIO_LIB && defined DISTRIBUTE
968 & brypiofile(idu3bc(iwest),ng), &
969# endif
970 & nbcfiles(ng), bry(1,ng), &
971 & recordless, update(1), &
972 & jlb, jub, n(ng), 2, 0, mm(ng)+1, n(ng), &
973 & boundary(ng) % uG_west)
974 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
975 END IF
976!
977 IF (tl_lbc(iwest,isvvel,ng)%acquire) THEN
978 CALL get_ngfld (ng, irpm, idv3bc(iwest), &
979 & bryncid(idv3bc(iwest),ng), &
980# if defined PIO_LIB && defined DISTRIBUTE
981 & brypiofile(idv3bc(iwest),ng), &
982# endif
983 & nbcfiles(ng), bry(1,ng), &
984 & recordless, update(1), &
985 & jlb, jub, n(ng), 2, 1, mm(ng)+1, n(ng), &
986 & boundary(ng) % vG_west)
987 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
988 END IF
989!
990 IF (tl_lbc(ieast,isuvel,ng)%acquire) THEN
991 CALL get_ngfld (ng, irpm, idu3bc(ieast), &
992 & bryncid(idu3bc(ieast),ng), &
993# if defined PIO_LIB && defined DISTRIBUTE
994 & brypiofile(idu3bc(ieast),ng), &
995# endif
996 & nbcfiles(ng), bry(1,ng), &
997 & recordless, update(1), &
998 & jlb, jub, n(ng), 2, 0, mm(ng)+1, n(ng), &
999 & boundary(ng) % uG_east)
1000 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1001 END IF
1002!
1003 IF (tl_lbc(ieast,isvvel,ng)%acquire) THEN
1004 CALL get_ngfld (ng, irpm, idv3bc(ieast), &
1005 & bryncid(idv3bc(ieast),ng), &
1006# if defined PIO_LIB && defined DISTRIBUTE
1007 & brypiofile(idv3bc(ieast),ng), &
1008# endif
1009 & nbcfiles(ng), bry(1,ng), &
1010 & recordless, update(1), &
1011 & jlb, jub, n(ng), 2, 1, mm(ng)+1, n(ng), &
1012 & boundary(ng) % vG_east)
1013 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1014 END IF
1015!
1016 IF (tl_lbc(isouth,isuvel,ng)%acquire) THEN
1017 CALL get_ngfld (ng, irpm, idu3bc(isouth), &
1018 & bryncid(idu3bc(isouth),ng), &
1019# if defined PIO_LIB && defined DISTRIBUTE
1020 & brypiofile(idu3bc(isouth),ng), &
1021# endif
1022 & nbcfiles(ng), bry(1,ng), &
1023 & recordless, update(1), &
1024 & ilb, iub, n(ng), 2, 1, lm(ng)+1, n(ng), &
1025 & boundary(ng) % uG_south)
1026 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1027 END IF
1028!
1029 IF (tl_lbc(isouth,isvvel,ng)%acquire) THEN
1030 CALL get_ngfld (ng, irpm, idv3bc(isouth), &
1031 & bryncid(idv3bc(isouth),ng), &
1032# if defined PIO_LIB && defined DISTRIBUTE
1033 & brypiofile(idv3bc(isouth),ng), &
1034# endif
1035 & nbcfiles(ng), bry(1,ng), &
1036 & recordless, update(1), &
1037 & ilb, iub, n(ng), 2, 0, lm(ng)+1, n(ng), &
1038 & boundary(ng) % vG_south)
1039 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1040 END IF
1041!
1042 IF (tl_lbc(inorth,isuvel,ng)%acquire) THEN
1043 CALL get_ngfld (ng, irpm, idu3bc(inorth), &
1044 & bryncid(idu3bc(inorth),ng), &
1045# if defined PIO_LIB && defined DISTRIBUTE
1046 & brypiofile(idu3bc(inorth),ng), &
1047# endif
1048 & nbcfiles(ng), bry(1,ng), &
1049 & recordless, update(1), &
1050 & ilb, iub, n(ng), 2, 1, lm(ng)+1, n(ng), &
1051 & boundary(ng) % uG_north)
1052 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1053 END IF
1054!
1055 IF (tl_lbc(inorth,isvvel,ng)%acquire) THEN
1056 CALL get_ngfld (ng, irpm, idv3bc(inorth), &
1057 & bryncid(idv3bc(inorth),ng), &
1058# if defined PIO_LIB && defined DISTRIBUTE
1059 & brypiofile(idv3bc(inorth),ng), &
1060# endif
1061 & nbcfiles(ng), bry(1,ng), &
1062 & recordless, update(1), &
1063 & ilb, iub, n(ng), 2, 0, lm(ng)+1, n(ng), &
1064 & boundary(ng) % vG_north)
1065 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1066 END IF
1067 END IF
1068# endif
1069
1070# ifndef ANA_TOBC
1071!
1072 IF (lprocessobc(ng)) THEN
1073 DO i=1,nt(ng)
1074 IF (tl_lbc(iwest,istvar(i),ng)%acquire) THEN
1075 CALL get_ngfld (ng, irpm, idtbry(iwest,i), &
1076 & bryncid(idtbry(iwest,i),ng), &
1077# if defined PIO_LIB && defined DISTRIBUTE
1078 & brypiofile(idtbry(iwest,i),ng), &
1079# endif
1080 & nbcfiles(ng), bry(1,ng), &
1081 & recordless, update(1), &
1082 & jlb, jub, n(ng), 2, 0, mm(ng)+1, n(ng), &
1083 & boundary(ng) % tG_west(:,:,:,i))
1084 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1085 END IF
1086 END DO
1087!
1088 DO i=1,nt(ng)
1089 IF (tl_lbc(ieast,istvar(i),ng)%acquire) THEN
1090 CALL get_ngfld (ng, irpm, idtbry(ieast,i), &
1091 & bryncid(idtbry(ieast,i),ng), &
1092# if defined PIO_LIB && defined DISTRIBUTE
1093 & brypiofile(idtbry(ieast,i),ng), &
1094# endif
1095 & nbcfiles(ng), bry(1,ng), &
1096 & recordless, update(1), &
1097 & jlb, jub, n(ng), 2, 0, mm(ng)+1, n(ng), &
1098 & boundary(ng) % tG_east(:,:,:,i))
1099 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1100 END IF
1101 END DO
1102!
1103 DO i=1,nt(ng)
1104 IF (tl_lbc(isouth,istvar(i),ng)%acquire) THEN
1105 CALL get_ngfld (ng, irpm, idtbry(isouth,i), &
1106 & bryncid(idtbry(isouth,i),ng), &
1107# if defined PIO_LIB && defined DISTRIBUTE
1108 & brypiofile(idtbry(isouth,i),ng), &
1109# endif
1110 & nbcfiles(ng), bry(1,ng), &
1111 & recordless, update(1), &
1112 & ilb, iub, n(ng), 2, 0, lm(ng)+1, n(ng), &
1113 & boundary(ng) % tG_south(:,:,:,i))
1114 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1115 END IF
1116 END DO
1117!
1118 DO i=1,nt(ng)
1119 IF (tl_lbc(inorth,istvar(i),ng)%acquire) THEN
1120 CALL get_ngfld (ng, irpm, idtbry(inorth,i), &
1121 & bryncid(idtbry(inorth,i),ng), &
1122# if defined PIO_LIB && defined DISTRIBUTE
1123 & brypiofile(idtbry(inorth,i),ng), &
1124# endif
1125 & nbcfiles(ng), bry(1,ng), &
1126 & recordless, update(1), &
1127 & ilb, iub, n(ng), 2, 0, lm(ng)+1, n(ng), &
1128 & boundary(ng) % tG_north(:,:,:,i))
1129 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1130 END IF
1131 END DO
1132 END IF
1133# endif
1134# endif
1135!
1136!=======================================================================
1137! Read in data from Climatology NetCDF file.
1138!=======================================================================
1139
1140# ifndef ANA_SSH
1141!
1142! Free-surface.
1143!
1144 IF (lsshclm(ng)) THEN
1145 CALL get_2dfld (ng, irpm, idsshc, clmncid(idsshc,ng), &
1146# if defined PIO_LIB && defined DISTRIBUTE
1147 & clmpiofile(idsshc,ng), &
1148# endif
1149 & nclmfiles(ng), clm(1,ng), update(1), &
1150 & lbi, ubi, lbj, ubj, 2, 1, &
1151# ifdef MASKING
1152 & grid(ng) % rmask, &
1153# endif
1154 & clima(ng) % sshG)
1155 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1156 END IF
1157# endif
1158# ifndef ANA_M2CLIMA
1159!
1160! 2D momentum.
1161!
1162 IF (lm2clm(ng)) THEN
1163 CALL get_2dfld (ng, irpm, idubcl, clmncid(idubcl,ng), &
1164# if defined PIO_LIB && defined DISTRIBUTE
1165 & clmpiofile(idubcl,ng), &
1166# endif
1167 & nclmfiles(ng), clm(1,ng), update(1), &
1168 & lbi, ubi, lbj, ubj, 2, 1, &
1169# ifdef MASKING
1170 & grid(ng) % umask, &
1171# endif
1172 & clima(ng) % ubarclmG)
1173 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1174!
1175 CALL get_2dfld (ng, irpm, idvbcl, clmncid(idvbcl,ng), &
1176# if defined PIO_LIB && defined DISTRIBUTE
1177 & clmpiofile(idvbcl,ng), &
1178# endif
1179 & nclmfiles(ng), clm(1,ng), update(1), &
1180 & lbi, ubi, lbj, ubj, 2, 1, &
1181# ifdef MASKING
1182 & grid(ng) % vmask, &
1183# endif
1184 & clima(ng) % vbarclmG)
1185 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1186 END IF
1187# endif
1188# ifdef SOLVE3D
1189# ifndef ANA_M3CLIMA
1190!
1191! 3D momentum.
1192!
1193 IF (lm3clm(ng)) THEN
1194 CALL get_3dfld (ng, irpm, iduclm, clmncid(iduclm,ng), &
1195# if defined PIO_LIB && defined DISTRIBUTE
1196 & clmpiofile(iduclm,ng), &
1197# endif
1198 & nclmfiles(ng), clm(1,ng), update(1), &
1199 & lbi, ubi, lbj, ubj, 1, n(ng), 2, 1, &
1200# ifdef MASKING
1201 & grid(ng) % umask, &
1202# endif
1203 & clima(ng) % uclmG)
1204 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1205!
1206 CALL get_3dfld (ng, irpm, idvclm, clmncid(idvclm,ng), &
1207# if defined PIO_LIB && defined DISTRIBUTE
1208 & clmpiofile(idvclm,ng), &
1209# endif
1210 & nclmfiles(ng), clm(1,ng), update(1), &
1211 & lbi, ubi, lbj, ubj, 1, n(ng), 2, 1, &
1212# ifdef MASKING
1213 & grid(ng) % vmask, &
1214# endif
1215 & clima(ng) % vclmG)
1216 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1217 END IF
1218# endif
1219# ifndef ANA_TCLIMA
1220!
1221! Tracers.
1222!
1223 ic=0
1224 DO i=1,nt(ng)
1225 IF (ltracerclm(i,ng)) THEN
1226 ic=ic+1
1227 CALL get_3dfld (ng, irpm, idtclm(i), &
1228 & clmncid(idtclm(i),ng), &
1229# if defined PIO_LIB && defined DISTRIBUTE
1230 & clmpiofile(idtclm(i),ng), &
1231# endif
1232 & nclmfiles(ng), clm(1,ng), update(1), &
1233 & lbi, ubi, lbj, ubj, 1, n(ng), 2, 1, &
1234# ifdef MASKING
1235 & grid(ng) % rmask, &
1236# endif
1237 & clima(ng) % tclmG(:,:,:,:,ic))
1238 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1239 END IF
1240 END DO
1241# endif
1242# endif
1243
1244# ifdef FORWARD_READ
1245!
1246!-----------------------------------------------------------------------
1247! Read in forward state solution.
1248!-----------------------------------------------------------------------
1249!
1250! Read in free-surface.
1251!
1252 CALL get_2dfld (ng, irpm, idfsur, fwd(ng)%ncid, &
1253# if defined PIO_LIB && defined DISTRIBUTE
1254 & fwd(ng)%pioFile, &
1255# endif
1256 & 1, fwd(ng), update(1), &
1257 & lbi, ubi, lbj, ubj, 2, 1, &
1258# ifdef MASKING
1259 & grid(ng) % rmask, &
1260# endif
1261 & ocean(ng) % zetaG)
1262 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1263!
1264! Read 2D momentum.
1265!
1266 CALL get_2dfld (ng, irpm, idubar, fwd(ng)%ncid, &
1267# if defined PIO_LIB && defined DISTRIBUTE
1268 & fwd(ng)%pioFile, &
1269# endif
1270 & 1, fwd(ng), update(1), &
1271 & lbi, ubi, lbj, ubj, 2, 1, &
1272# ifdef MASKING
1273 & grid(ng) % umask, &
1274# endif
1275 & ocean(ng) % ubarG)
1276 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1277
1278 CALL get_2dfld (ng, irpm, idvbar, fwd(ng)%ncid, &
1279# if defined PIO_LIB && defined DISTRIBUTE
1280 & fwd(ng)%pioFile, &
1281# endif
1282 & 1, fwd(ng), update(1), &
1283 & lbi, ubi, lbj, ubj, 2, 1, &
1284# ifdef MASKING
1285 & grid(ng) % vmask, &
1286# endif
1287 & ocean(ng) % vbarG)
1288 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1289
1290# ifdef FORWARD_RHS
1291!
1292! Read in variables associated with 2D right-hand-side terms.
1293!
1294 CALL get_2dfld (ng, irpm, idrzet, fwd(ng)%ncid, &
1295# if defined PIO_LIB && defined DISTRIBUTE
1296 & fwd(ng)%pioFile, &
1297# endif
1298 & 1, fwd(ng), update(1), &
1299 & lbi, ubi, lbj, ubj, 2, 1, &
1300# ifdef MASKING
1301 & grid(ng) % rmask, &
1302# endif
1303 & ocean(ng) % rzetaG)
1304 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1305!
1306 CALL get_2dfld (ng, irpm, idru2d, fwd(ng)%ncid, &
1307# if defined PIO_LIB && defined DISTRIBUTE
1308 & fwd(ng)%pioFile, &
1309# endif
1310 & 1, fwd(ng), update(1), &
1311 & lbi, ubi, lbj, ubj, 2, 1, &
1312# ifdef MASKING
1313 & grid(ng) % umask, &
1314# endif
1315 & ocean(ng) % rubarG)
1316 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1317!
1318 CALL get_2dfld (ng, irpm, idrv2d, fwd(ng)%ncid, &
1319# if defined PIO_LIB && defined DISTRIBUTE
1320 & fwd(ng)%pioFile, &
1321# endif
1322 & 1, fwd(ng), update(1), &
1323 & lbi, ubi, lbj, ubj, 2, 1, &
1324# ifdef MASKING
1325 & grid(ng) % vmask, &
1326# endif
1327 & ocean(ng) % rvbarG)
1328 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1329# endif
1330
1331# ifdef SOLVE3D
1332!
1333! Read in variables associated with time-averaged 2D momentum terms.
1334!
1335 CALL get_2dfld (ng, irpm, idufx1, fwd(ng)%ncid, &
1336# if defined PIO_LIB && defined DISTRIBUTE
1337 & fwd(ng)%pioFile, &
1338# endif
1339 & 1, fwd(ng), update(1), &
1340 & lbi, ubi, lbj, ubj, 2, 1, &
1341# ifdef MASKING
1342 & grid(ng) % umask, &
1343# endif
1344 & coupling(ng) % DU_avg1G)
1345 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1346!
1347 CALL get_2dfld (ng, irpm, idufx2, fwd(ng)%ncid, &
1348# if defined PIO_LIB && defined DISTRIBUTE
1349 & fwd(ng)%pioFile, &
1350# endif
1351 & 1, fwd(ng), update(1), &
1352 & lbi, ubi, lbj, ubj, 2, 1, &
1353# ifdef MASKING
1354 & grid(ng) % umask, &
1355# endif
1356 & coupling(ng) % DU_avg2G)
1357 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1358!
1359 CALL get_2dfld (ng, irpm, idvfx1, fwd(ng)%ncid, &
1360# if defined PIO_LIB && defined DISTRIBUTE
1361 & fwd(ng)%pioFile, &
1362# endif
1363 & 1, fwd(ng), update(1), &
1364 & lbi, ubi, lbj, ubj, 2, 1, &
1365# ifdef MASKING
1366 & grid(ng) % vmask, &
1367# endif
1368 & coupling(ng) % DV_avg1G)
1369 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1370!
1371 CALL get_2dfld (ng, irpm, idvfx2, fwd(ng)%ncid, &
1372# if defined PIO_LIB && defined DISTRIBUTE
1373 & fwd(ng)%pioFile, &
1374# endif
1375 & 1, fwd(ng), update(1), &
1376 & lbi, ubi, lbj, ubj, 2, 1, &
1377# ifdef MASKING
1378 & grid(ng) % vmask, &
1379# endif
1380 & coupling(ng) % DV_avg2G)
1381 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1382!
1383! Read in 3D momentum.
1384!
1385 CALL get_3dfld (ng, irpm, iduvel, fwd(ng)%ncid, &
1386# if defined PIO_LIB && defined DISTRIBUTE
1387 & fwd(ng)%pioFile, &
1388# endif
1389 & 1, fwd(ng), update(1), &
1390 & lbi, ubi, lbj, ubj, 1, n(ng), 2, 1, &
1391# ifdef MASKING
1392 & grid(ng) % umask, &
1393# endif
1394 & ocean(ng) % uG)
1395 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1396!
1397 CALL get_3dfld (ng, irpm, idvvel, fwd(ng)%ncid, &
1398# if defined PIO_LIB && defined DISTRIBUTE
1399 & fwd(ng)%pioFile, &
1400# endif
1401 & 1, fwd(ng), update(1), &
1402 & lbi, ubi, lbj, ubj, 1, n(ng), 2, 1, &
1403# ifdef MASKING
1404 & grid(ng) % vmask, &
1405# endif
1406 & ocean(ng) % vG)
1407 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1408
1409# ifdef FORWARD_RHS
1410!
1411! Read in variables associated with 3D momentum right-hand-side terms.
1412!
1413 CALL get_2dfld (ng, irpm, idruct, fwd(ng)%ncid, &
1414# if defined PIO_LIB && defined DISTRIBUTE
1415 & fwd(ng)%pioFile, &
1416# endif
1417 & 1, fwd(ng), update(1), &
1418 & lbi, ubi, lbj, ubj, 2, 1, &
1419# ifdef MASKING
1420 & grid(ng) % umask, &
1421# endif
1422 & coupling(ng) % rufrcG)
1423 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1424!
1425 CALL get_2dfld (ng, irpm, idrvct, fwd(ng)%ncid, &
1426# if defined PIO_LIB && defined DISTRIBUTE
1427 & fwd(ng)%pioFile, &
1428# endif
1429 & 1, fwd(ng), update(1), &
1430 & lbi, ubi, lbj, ubj, 2, 1, &
1431# ifdef MASKING
1432 & grid(ng) % vmask, &
1433# endif
1434 & coupling(ng) % rvfrcG)
1435 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1436!
1437 CALL get_3dfld (ng, irpm, idru3d, fwd(ng)%ncid, &
1438# if defined PIO_LIB && defined DISTRIBUTE
1439 & fwd(ng)%pioFile, &
1440# endif
1441 & 1, fwd(ng), update(1), &
1442 & lbi, ubi, lbj, ubj, 1, n(ng), 2, 1, &
1443# ifdef MASKING
1444 & grid(ng) % umask, &
1445# endif
1446 & ocean(ng) % ruG)
1447 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1448!
1449 CALL get_3dfld (ng, irpm, idrv3d, fwd(ng)%ncid, &
1450# if defined PIO_LIB && defined DISTRIBUTE
1451 & fwd(ng)%pioFile, &
1452# endif
1453 & 1, fwd(ng), update(1), &
1454 & lbi, ubi, lbj, ubj, 1, n(ng), 2, 1, &
1455# ifdef MASKING
1456 & grid(ng) % vmask, &
1457# endif
1458 & ocean(ng) % rvG)
1459 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1460# endif
1461!
1462! Read in 3D tracers.
1463!
1464 DO i=1,nt(ng)
1465 CALL get_3dfld (ng, irpm, idtvar(i), fwd(ng)%ncid, &
1466# if defined PIO_LIB && defined DISTRIBUTE
1467 & fwd(ng)%pioFile, &
1468# endif
1469 & 1, fwd(ng), update(1), &
1470 & lbi, ubi, lbj, ubj, 1, n(ng), 2, 1, &
1471# ifdef MASKING
1472 & grid(ng) % rmask, &
1473# endif
1474 & ocean(ng) % tG(:,:,:,:,i))
1475 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1476 END DO
1477
1478# ifdef FORWARD_MIXING
1479!
1480! Read in vertical mixing variables.
1481!
1482 DO i=1,nat
1483 scale=fscale(iddiff(i),ng) ! save and rescale
1484 fscale(iddiff(i),ng)=tl_akt_fac(i,ng)
1485 CALL get_3dfld (ng, irpm, iddiff(i), fwd(ng)%ncid, &
1486# if defined PIO_LIB && defined DISTRIBUTE
1487 & fwd(ng)%pioFile, &
1488# endif
1489 & 1, fwd(ng), update(1), &
1490 & lbi, ubi, lbj, ubj, 0, n(ng), 2, 1, &
1491# ifdef MASKING
1492 & grid(ng) % rmask, &
1493# endif
1494 & mixing(ng) % AktG(:,:,:,:,i))
1495 fscale(iddiff(i),ng)=scale
1496 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1497 END DO
1498!
1499 scale=fscale(idvvis,ng) ! save and rescale
1500 fscale(idvvis,ng)=tl_akv_fac(ng)
1501 CALL get_3dfld (ng, irpm, idvvis, fwd(ng)%ncid, &
1502# if defined PIO_LIB && defined DISTRIBUTE
1503 & fwd(ng)%pioFile, &
1504# endif
1505 & 1, fwd(ng), update(1), &
1506 & lbi, ubi, lbj, ubj, 0, n(ng), 2, 1, &
1507# ifdef MASKING
1508 & grid(ng) % rmask, &
1509# endif
1510 & mixing(ng) % AkvG)
1511 fscale(idvvis,ng)=scale
1512 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1513# endif
1514
1515# if defined MY25_MIXING_NOT_YET || defined GLS_MIXING_NOT_YET
1516!
1517! Read in turbulent kinetic energy.
1518!
1519 CALL get_3dfld (ng, irpm, idmtke, fwd(ng)%ncid, &
1520# if defined PIO_LIB && defined DISTRIBUTE
1521 & fwd(ng)%pioFile, &
1522# endif
1523 & 1, fwd(ng), update(1), &
1524 & lbi, ubi, lbj, ubj, 0, n(ng), 2, 1, &
1525# ifdef MASKING
1526 & grid(ng) % rmask, &
1527# endif
1528 & mixing(ng) % tkeG)
1529 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1530!
1531! Read in turbulent kinetic energy times length scale.
1532!
1533 CALL get_3dfld (ng, irpm, idmtls, fwd(ng)%ncid, &
1534# if defined PIO_LIB && defined DISTRIBUTE
1535 & fwd(ng)%pioFile, &
1536# endif
1537 & 1, fwd(ng), update(1), &
1538 & lbi, ubi, lbj, ubj, 0, n(ng), 2, 1, &
1539# ifdef MASKING
1540 & grid(ng) % rmask, &
1541# endif
1542 & mixing(ng) % glsG)
1543 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1544!
1545! Read in vertical mixing length scale.
1546!
1547 CALL get_3dfld (ng, irpm, idvmls, fwd(ng)%ncid, &
1548# if defined PIO_LIB && defined DISTRIBUTE
1549 & fwd(ng)%pioFile, &
1550# endif
1551 & 1, fwd(ng), update(1), &
1552 & lbi, ubi, lbj, ubj, 0, n(ng), 2, 1, &
1553# ifdef MASKING
1554 & grid(ng) % rmask, &
1555# endif
1556 & mixing(ng) % LscaleG)
1557 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1558!
1559! Read in vertical mixing coefficient for turbulent kinetic energy.
1560!
1561 CALL get_3dfld (ng, irpm, idvmkk, fwd(ng)%ncid, &
1562# if defined PIO_LIB && defined DISTRIBUTE
1563 & fwd(ng)%pioFile, &
1564# endif
1565 & 1, fwd(ng), update(1), &
1566 & lbi, ubi, lbj, ubj, 0, n(ng), 2, 1, &
1567# ifdef MASKING
1568 & grid(ng) % rmask, &
1569# endif
1570 & mixing(ng) % AkkG)
1571 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1572
1573# ifdef GLS_MIXING_NOT_YET
1574!
1575! Read in vertical mixing coefficient for turbulent length scale.
1576!
1577 CALL get_3dfld (ng, irpm, idvmkp, fwd(ng)%ncid, &
1578# if defined PIO_LIB && defined DISTRIBUTE
1579 & fwd(ng)%pioFile, &
1580# endif
1581 & 1, fwd(ng), update(1), &
1582 & lbi, ubi, lbj, ubj, 0, n(ng), 2, 1, &
1583# ifdef MASKING
1584 & grid(ng) % rmask, &
1585# endif
1586 & mixing(ng) % AkpG)
1587 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1588# endif
1589# endif
1590
1591# ifdef LMD_MIXING_NOT_YET
1592!
1593! Read in depth of surface oceanic boundary layer.
1594!
1595 CALL get_2dfld (ng, irpm, idhsbl, fwd(ng)%ncid, &
1596# if defined PIO_LIB && defined DISTRIBUTE
1597 & fwd(ng)%pioFile, &
1598# endif
1599 & 1, fwd(ng), update(1), &
1600 & lbi, ubi, lbj, ubj, 2, 1, &
1601# ifdef MASKING
1602 & grid(ng) % rmask, &
1603# endif
1604 & mixing(ng) % hsblG)
1605 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1606# endif
1607
1608# ifdef LMD_BKPP_NOT_YET
1609!
1610! Read in depth of bottom oceanic boundary layer.
1611!
1612 CALL get_2dfld (ng, irpm, idhbbl, fwd(ng)%ncid, &
1613# if defined PIO_LIB && defined DISTRIBUTE
1614 & fwd(ng)%pioFile, &
1615# endif
1616 & 1, fwd(ng), update(1), &
1617 & lbi, ubi, lbj, ubj, 2, 1, &
1618# ifdef MASKING
1619 & grid(ng) % rmask, &
1620# endif
1621 & mixing(ng) % hbblG)
1622 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1623# endif
1624
1625# ifdef LMD_NONLOCAL_NOT_YET
1626!
1627! Read in boundary layer nonlocal transport.
1628!
1629 DO i=1,nat
1630 CALL get_3dfld (ng, irpm, idghat(i), fwd(ng)%ncid, &
1631# if defined PIO_LIB && defined DISTRIBUTE
1632 & fwd(ng)%pioFile, &
1633# endif
1634 & 1, fwd(ng), update(1), &
1635 & lbi, ubi, lbj, ubj, 0, n(ng), 2, 1, &
1636# ifdef MASKING
1637 & grid(ng) % rmask, &
1638# endif
1639 & mixing(ng) % ghatsG(:,:,:,i))
1640 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1641 END DO
1642# endif
1643# endif
1644
1645# if defined R4DVAR || defined R4DVAR_ANA_SENSITIVITY || \
1646 defined tl_r4dvar
1647!
1648!-----------------------------------------------------------------------
1649! Read frequent impulse forcing for weak constraint.
1650!-----------------------------------------------------------------------
1651!
1652 IF (frequentimpulse(ng)) THEN
1653 CALL get_2dfld (ng, irpm, 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, irpm, 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, irpm, 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, irpm, 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, irpm, 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, irpm, 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, irpm, 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 rp_initial(), and rp_main3d().

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