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