28# if defined PIO_LIB && defined DISTRIBUTE
43 PRIVATE :: get_extract_nf90
44# if defined PIO_LIB && defined DISTRIBUTE
45 PRIVATE :: get_extract_pio
51 SUBROUTINE get_extract (ng, tile, model)
56 integer,
intent(in) :: ng, tile, model
60 integer :: LBi, UBi, LBj, UBj
62 character (len=*),
parameter :: MyFile = &
69 lbi=xtr_bounds(ng)%LBi(tile)
70 ubi=xtr_bounds(ng)%UBi(tile)
71 lbj=xtr_bounds(ng)%LBj(tile)
72 ubj=xtr_bounds(ng)%UBj(tile)
74 SELECT CASE (
grx(ng)%IOtype)
76 CALL get_extract_nf90 (ng, tile, model, &
79# if defined PIO_LIB && defined DISTRIBUTE
81 CALL get_extract_pio (ng, tile, model, &
88 IF (
founderror(exit_flag, noerror, __line__, myfile))
RETURN
95 IF (extractflag(ng).ge.1)
THEN
96 CALL interp_coords (ng, tile, model,
p2dvar, &
97 &
grid(ng)%Gangle_psi, &
99 &
grid(ng)%Gmask_psi, &
104 & extract(ng)%Gmask_psi, &
106 & extract(ng)%Gx_psi, &
107 & extract(ng)%Gy_psi, &
108 & extract(ng)%Iout_psi, &
109 & extract(ng)%Jout_psi)
111 CALL interp_coords (ng, tile, model,
r2dvar, &
112 &
grid(ng)%Gangle_rho, &
114 &
grid(ng)%Gmask_rho, &
119 & extract(ng)%Gmask_rho, &
121 & extract(ng)%Gx_rho, &
122 & extract(ng)%Gy_rho, &
123 & extract(ng)%Iout_rho, &
124 & extract(ng)%Jout_rho)
126 CALL interp_coords (ng, tile, model,
u2dvar, &
127 &
grid(ng)%Gangle_u, &
129 &
grid(ng)%Gmask_u, &
134 & extract(ng)%Gmask_u, &
136 & extract(ng)%Gx_u, &
137 & extract(ng)%Gy_u, &
138 & extract(ng)%Iout_u, &
139 & extract(ng)%Jout_u)
141 CALL interp_coords (ng, tile, model,
v2dvar, &
142 &
grid(ng)%Gangle_v, &
144 &
grid(ng)%Gmask_v, &
149 & extract(ng)%Gmask_v, &
151 & extract(ng)%Gx_v, &
152 & extract(ng)%Gy_v, &
154 & extract(ng)%Iout_v, &
155 & extract(ng)%Jout_v)
158 10
FORMAT (
' GET_EXTRACT - Illegal input file type, io_type = ',i0, &
159 & /,12x,
'Check KeyWord ''INP_LIB'' in ''roms.in''.')
162 END SUBROUTINE get_extract
165 SUBROUTINE get_extract_nf90 (ng, tile, model, &
166 & LBi, UBi, LBj, UBj)
171 integer,
intent(in) :: ng, tile, model
172 integer,
intent(in) :: LBi, UBi, LBj, UBj
176 integer :: cr, gtype, i, status, vindex
179 integer(i8b) :: Fhash
182 real(dp),
parameter :: Fscl = 1.0_dp
184 real(r8) :: Fmax, Fmin
186 character (len=256) :: ncname
188 character (len=*),
parameter :: MyFile = &
189 & __FILE__//
", get_extract_nf90"
198 IF (
founderror(exit_flag, noerror, __line__, myfile))
RETURN
203 IF (
grx(ng)%ncid.eq.-1)
THEN
205 IF (
founderror(exit_flag, noerror, __line__, myfile))
THEN
206 WRITE (
stdout,10) trim(ncname)
214 & ncid =
grx(ng)%ncid)
215 IF (
founderror(exit_flag, noerror, __line__, myfile))
RETURN
241# if (defined CURVGRID && defined UV_ADV)
289 & ncid =
grx(ng)%ncid)
290 IF (
founderror(exit_flag, noerror, __line__, myfile))
RETURN
310 SELECT CASE (trim(adjustl(
var_name(i))))
316 status=nf_fread2d_xtr(ng, model, ncname,
grx(ng)%ncid, &
319 & lbi, ubi, lbj, ubj, &
320 & fscl, fmin, fmax, &
322 & extract(ng) % rmask, &
330 IF (
founderror(status, nf90_noerr, __line__, myfile))
THEN
336 WRITE (
stdout,30)
'bathymetry at RHO-points: h', &
337 & ng, trim(ncname), fmin, fmax
343 IF (ewperiodic(ng).or.nsperiodic(ng))
THEN
344 CALL exchange_r2d_xtr_tile (ng, tile, &
345 & lbi, ubi, lbj, ubj, &
349 CALL mp_exchange2d_xtr (ng, tile, model, 1, &
350 & lbi, ubi, lbj, ubj, &
352 & ewperiodic(ng), nsperiodic(ng), &
361 status=nf_fread2d_xtr(ng, model, ncname,
grx(ng)%ncid, &
364 & lbi, ubi, lbj, ubj, &
365 & fscl, fmin, fmax, &
366 & extract(ng) % rmask, &
368 & extract(ng) % rmask, &
371 & extract(ng) % rmask)
373 IF (
founderror(status, nf90_noerr, __line__, myfile))
THEN
379 WRITE (
stdout,30)
'mask on RHO-points: mask_rho', &
380 & ng, trim(ncname), fmin, fmax
386 IF (ewperiodic(ng).or.nsperiodic(ng))
THEN
387 CALL exchange_r2d_xtr_tile (ng, tile, &
388 & lbi, ubi, lbj, ubj, &
389 & extract(ng) % rmask)
392 CALL mp_exchange2d_xtr (ng, tile, model, 1, &
393 & lbi, ubi, lbj, ubj, &
395 & ewperiodic(ng), nsperiodic(ng), &
396 & extract(ng) % rmask)
403 status=nf_fread2d_xtr(ng, model, ncname,
grx(ng)%ncid, &
406 & lbi, ubi, lbj, ubj, &
407 & fscl, fmin, fmax, &
408 & extract(ng) % umask, &
410 & extract(ng) % umask, &
413 & extract(ng) % umask)
415 IF (
founderror(status, nf90_noerr, __line__, myfile))
THEN
421 WRITE (
stdout,30)
'mask on U-points: mask_u', &
422 & ng, trim(ncname), fmin, fmax
428 IF (ewperiodic(ng).or.nsperiodic(ng))
THEN
429 CALL exchange_u2d_xtr_tile (ng, tile, &
430 & lbi, ubi, lbj, ubj, &
431 & extract(ng) % umask)
434 CALL mp_exchange2d_xtr (ng, tile, model, 1, &
435 & lbi, ubi, lbj, ubj, &
437 & ewperiodic(ng), nsperiodic(ng), &
438 & extract(ng) % umask)
445 status=nf_fread2d_xtr(ng, model, ncname,
grx(ng)%ncid, &
448 & lbi, ubi, lbj, ubj, &
449 & fscl, fmin, fmax, &
450 & extract(ng) % vmask, &
452 & extract(ng) % vmask, &
455 & extract(ng) % vmask)
457 IF (
founderror(status, nf90_noerr, __line__, myfile))
THEN
463 WRITE (
stdout,30)
'mask on V-points: mask_v', &
464 & ng, trim(ncname), fmin, fmax
470 IF (ewperiodic(ng).or.nsperiodic(ng))
THEN
471 CALL exchange_v2d_xtr_tile (ng, tile, &
472 & lbi, ubi, lbj, ubj, &
473 & extract(ng) % vmask)
476 CALL mp_exchange2d_xtr (ng, tile, model, 1, &
477 & lbi, ubi, lbj, ubj, &
479 & ewperiodic(ng), nsperiodic(ng), &
480 & extract(ng) % vmask)
487 status=nf_fread2d_xtr(ng, model, ncname,
grx(ng)%ncid, &
490 & lbi, ubi, lbj, ubj, &
491 & fscl, fmin, fmax, &
492 & extract(ng) % pmask, &
494 & extract(ng) % pmask, &
497 & extract(ng) % pmask)
499 IF (
founderror(status, nf90_noerr, __line__, myfile))
THEN
505 WRITE (
stdout,30)
'mask on PSI-points: mask_psi', &
506 & ng, trim(ncname), fmin, fmax
512 IF (ewperiodic(ng).or.nsperiodic(ng))
THEN
513 CALL exchange_p2d_xtr_tile (ng, tile, &
514 & lbi, ubi, lbj, ubj, &
515 & extract(ng) % pmask)
518 CALL mp_exchange2d_xtr (ng, tile, model, 1, &
519 & lbi, ubi, lbj, ubj, &
521 & ewperiodic(ng), nsperiodic(ng), &
522 & extract(ng) % pmask)
530 status=nf_fread2d_xtr(ng, model, ncname,
grx(ng)%ncid, &
533 & lbi, ubi, lbj, ubj, &
534 & fscl, fmin, fmax, &
536 & extract(ng) % rmask, &
544 IF (
founderror(status, nf90_noerr, __line__, myfile))
THEN
550 WRITE (
stdout,30)
'Coriolis parameter at RHO-points: f',&
551 & ng, trim(ncname), fmin, fmax
557 IF (ewperiodic(ng).or.nsperiodic(ng))
THEN
558 CALL exchange_r2d_xtr_tile (ng, tile, &
559 & lbi, ubi, lbj, ubj, &
563 CALL mp_exchange2d_xtr (ng, tile, model, 1, &
564 & lbi, ubi, lbj, ubj, &
566 & ewperiodic(ng), nsperiodic(ng), &
575 status=nf_fread2d_xtr(ng, model, ncname,
grx(ng)%ncid, &
578 & lbi, ubi, lbj, ubj, &
579 & fscl, fmin, fmax, &
581 & extract(ng) % rmask, &
584 & extract(ng) % pm, &
589 IF (
founderror(status, nf90_noerr, __line__, myfile))
THEN
595 WRITE (
stdout,30)
'reciprocal XI-grid spacing: pm', &
596 & ng, trim(ncname), fmin, fmax
602 IF (ewperiodic(ng).or.nsperiodic(ng))
THEN
603 CALL exchange_r2d_xtr_tile (ng, tile, &
604 & lbi, ubi, lbj, ubj, &
608 CALL mp_exchange2d_xtr (ng, tile, model, 1, &
609 & lbi, ubi, lbj, ubj, &
611 & ewperiodic(ng), nsperiodic(ng), &
620 status=nf_fread2d_xtr(ng, model, ncname,
grx(ng)%ncid, &
623 & lbi, ubi, lbj, ubj, &
624 & fscl, fmin, fmax, &
626 & extract(ng) % rmask, &
629 & extract(ng) % pn, &
634 IF (
founderror(status, nf90_noerr, __line__, myfile))
THEN
640 WRITE (
stdout,30)
'reciprocal ETA-grid spacing: pn', &
641 & ng, trim(ncname), fmin, fmax
647 IF (ewperiodic(ng).or.nsperiodic(ng))
THEN
648 CALL exchange_r2d_xtr_tile (ng, tile, &
649 & lbi, ubi, lbj, ubj, &
653 CALL mp_exchange2d_xtr (ng, tile, model, 1, &
654 & lbi, ubi, lbj, ubj, &
656 & ewperiodic(ng), nsperiodic(ng), &
659# if (defined CURVGRID && defined UV_ADV)
665 status=nf_fread2d_xtr(ng, model, ncname,
grx(ng)%ncid, &
668 & lbi, ubi, lbj, ubj, &
669 & fscl, fmin, fmax, &
671 & extract(ng) % rmask, &
674 & extract(ng) % dmde, &
677 & extract(ng) % dmde)
679 IF (
founderror(status, nf90_noerr, __line__, myfile))
THEN
685 WRITE (
stdout,30)
'ETA-derivative of inverse metric '// &
686 &
'factor pm: dmde', &
687 & ng, trim(ncname), fmin, fmax
693 IF (ewperiodic(ng).or.nsperiodic(ng))
THEN
694 CALL exchange_r2d_xtr_tile (ng, tile, &
695 & lbi, ubi, lbj, ubj, &
696 & extract(ng) % dmde)
699 CALL mp_exchange2d_xtr (ng, tile, model, 1, &
700 & lbi, ubi, lbj, ubj, &
702 & ewperiodic(ng), nsperiodic(ng), &
703 & extract(ng) % dmde)
710 status=nf_fread2d_xtr(ng, model, ncname,
grx(ng)%ncid, &
713 & lbi, ubi, lbj, ubj, &
714 & fscl, fmin, fmax, &
716 & extract(ng) % rmask, &
719 & extract(ng) % dndx, &
722 & extract(ng) % dndx)
724 IF (
founderror(status, nf90_noerr, __line__, myfile))
THEN
730 WRITE (
stdout,30)
'XI-derivative of inverse metric '// &
731 &
'factor pn: dndx', &
732 & ng, trim(ncname), fmin, fmax
738 IF (ewperiodic(ng).or.nsperiodic(ng))
THEN
739 CALL exchange_r2d_xtr_tile (ng, tile, &
740 & lbi, ubi, lbj, ubj, &
741 & extract(ng) % dndx)
744 CALL mp_exchange2d_xtr (ng, tile, model, 1, &
745 & lbi, ubi, lbj, ubj, &
747 & ewperiodic(ng), nsperiodic(ng), &
748 & extract(ng) % dndx)
756 status=nf_fread2d_xtr(ng, model, ncname,
grx(ng)%ncid, &
759 & lbi, ubi, lbj, ubj, &
760 & fscl, fmin, fmax, &
762 & extract(ng) % pmask, &
765 & extract(ng) % xp, &
770 IF (
founderror(status, nf90_noerr, __line__, myfile))
THEN
776 WRITE (
stdout,30)
'x-location of PSI-points: x_psi', &
777 & ng, trim(ncname), fmin, fmax
784 CALL mp_exchange2d_xtr (ng, tile, model, 1, &
785 & lbi, ubi, lbj, ubj, &
787 & .false., .false., &
795 status=nf_fread2d_xtr(ng, model, ncname,
grx(ng)%ncid, &
798 & lbi, ubi, lbj, ubj, &
799 & fscl, fmin, fmax, &
801 & extract(ng) % pmask, &
804 & extract(ng) % yp, &
809 IF (
founderror(status, nf90_noerr, __line__, myfile))
THEN
815 WRITE (
stdout,30)
'y-location of PSI-points: y-psi', &
816 & ng, trim(ncname), fmin, fmax
823 CALL mp_exchange2d_xtr (ng, tile, model, 1, &
824 & lbi, ubi, lbj, ubj, &
826 & .false., .false., &
834 status=nf_fread2d_xtr(ng, model, ncname,
grx(ng)%ncid, &
837 & lbi, ubi, lbj, ubj, &
838 & fscl, fmin, fmax, &
840 & extract(ng) % rmask, &
843 & extract(ng) % xr, &
848 IF (
founderror(status, nf90_noerr, __line__, myfile))
THEN
854 WRITE (
stdout,30)
'x-location of RHO-points: x-rho', &
855 & ng, trim(ncname), fmin, fmax
861 IF (.not.spherical)
THEN
862 extract(ng)%LonMin(ng)=fmin
863 extract(ng)%LonMax(ng)=fmax
866 CALL mp_exchange2d_xtr (ng, tile, model, 1, &
867 & lbi, ubi, lbj, ubj, &
869 & .false., .false., &
877 status=nf_fread2d_xtr(ng, model, ncname,
grx(ng)%ncid, &
880 & lbi, ubi, lbj, ubj, &
881 & fscl, fmin, fmax, &
883 & extract(ng) % rmask, &
886 & extract(ng) % yr, &
891 IF (
founderror(status, nf90_noerr, __line__, myfile))
THEN
897 WRITE (
stdout,30)
'y-location of RHO-points: y_rho', &
898 & ng, trim(ncname), fmin, fmax
904 IF (.not.spherical)
THEN
905 extract(ng)%LatMin(ng)=fmin
906 extract(ng)%LatMax(ng)=fmax
909 CALL mp_exchange2d_xtr (ng, tile, model, 1, &
910 & lbi, ubi, lbj, ubj, &
912 & .false., .false., &
920 status=nf_fread2d_xtr(ng, model, ncname,
grx(ng)%ncid, &
923 & lbi, ubi, lbj, ubj, &
924 & fscl, fmin, fmax, &
926 & extract(ng) % umask, &
929 & extract(ng) % xu, &
934 IF (
founderror(status, nf90_noerr, __line__, myfile))
THEN
940 WRITE (
stdout,30)
'x-location of U-points: x_u', &
941 & ng, trim(ncname), fmin, fmax
948 CALL mp_exchange2d_xtr (ng, tile, model, 1, &
949 & lbi, ubi, lbj, ubj, &
951 & .false., .false., &
959 status=nf_fread2d_xtr(ng, model, ncname,
grx(ng)%ncid, &
962 & lbi, ubi, lbj, ubj, &
963 & fscl, fmin, fmax, &
965 & extract(ng) % umask, &
968 & extract(ng) % yu, &
973 IF (
founderror(status, nf90_noerr, __line__, myfile))
THEN
979 WRITE (
stdout,30)
'y-location of U-points: y_u', &
980 & ng, trim(ncname), fmin, fmax
987 CALL mp_exchange2d_xtr (ng, tile, model, 1, &
988 & lbi, ubi, lbj, ubj, &
990 & .false., .false., &
998 status=nf_fread2d_xtr(ng, model, ncname,
grx(ng)%ncid, &
1000 & 0, gtype, vsize, &
1001 & lbi, ubi, lbj, ubj, &
1002 & fscl, fmin, fmax, &
1004 & extract(ng) % vmask, &
1007 & extract(ng) % xv, &
1012 IF (
founderror(status, nf90_noerr, __line__, myfile))
THEN
1018 WRITE (
stdout,30)
'x-location of V-points: x_v', &
1019 & ng, trim(ncname), fmin, fmax
1026 CALL mp_exchange2d_xtr (ng, tile, model, 1, &
1027 & lbi, ubi, lbj, ubj, &
1029 & .false., .false., &
1037 status=nf_fread2d_xtr(ng, model, ncname,
grx(ng)%ncid, &
1039 & 0, gtype, vsize, &
1040 & lbi, ubi, lbj, ubj, &
1041 & fscl, fmin, fmax, &
1043 & extract(ng) % vmask, &
1046 & extract(ng) % yv, &
1051 IF (
founderror(status, nf90_noerr, __line__, myfile))
THEN
1057 WRITE (
stdout,30)
'y-location of V-points: y_v', &
1058 & ng, trim(ncname), fmin, fmax
1065 CALL mp_exchange2d_xtr (ng, tile, model, 1, &
1066 & lbi, ubi, lbj, ubj, &
1068 & .false., .false., &
1077 status=nf_fread2d_xtr(ng, model, ncname,
grx(ng)%ncid, &
1079 & 0, gtype, vsize, &
1080 & lbi, ubi, lbj, ubj, &
1081 & fscl, fmin, fmax, &
1083 & extract(ng) % pmask, &
1086 & extract(ng) % lonp, &
1089 & extract(ng) % lonp)
1091 IF (
founderror(status, nf90_noerr, __line__, myfile))
THEN
1097 WRITE (
stdout,30)
'longitude of PSI-points: lon_psi', &
1098 & ng, trim(ncname), fmin, fmax
1105 CALL mp_exchange2d_xtr (ng, tile, model, 1, &
1106 & lbi, ubi, lbj, ubj, &
1108 & .false., .false., &
1109 & extract(ng) % lonp)
1118 status=nf_fread2d_xtr(ng, model, ncname,
grx(ng)%ncid, &
1120 & 0, gtype, vsize, &
1121 & lbi, ubi, lbj, ubj, &
1122 & fscl, fmin, fmax, &
1124 & extract(ng) % pmask, &
1127 & extract(ng) % latp, &
1130 & extract(ng) % latp)
1132 IF (
founderror(status, nf90_noerr, __line__, myfile))
THEN
1138 WRITE (
stdout,30)
'latitude of PSI-points lat_psi', &
1139 & ng, trim(ncname), fmin, fmax
1146 CALL mp_exchange2d_xtr (ng, tile, model, 1, &
1147 & lbi, ubi, lbj, ubj, &
1149 & .false., .false., &
1150 & extract(ng) % latp)
1159 status=nf_fread2d_xtr(ng, model, ncname,
grx(ng)%ncid, &
1161 & 0, gtype, vsize, &
1162 & lbi, ubi, lbj, ubj, &
1163 & fscl, fmin, fmax, &
1165 & extract(ng) % rmask, &
1168 & extract(ng) % lonr, &
1171 & extract(ng) % lonr)
1173 IF (
founderror(status, nf90_noerr, __line__, myfile))
THEN
1179 WRITE (
stdout,30)
'longitude of RHO-points: lon_rho', &
1180 & ng, trim(ncname), fmin, fmax
1186 extract(ng)%LonMin(ng)=fmin
1187 extract(ng)%LonMax(ng)=fmax
1190 CALL mp_exchange2d_xtr (ng, tile, model, 1, &
1191 & lbi, ubi, lbj, ubj, &
1193 & .false., .false., &
1194 & extract(ng) % lonr)
1203 status=nf_fread2d_xtr(ng, model, ncname,
grx(ng)%ncid, &
1205 & 0, gtype, vsize, &
1206 & lbi, ubi, lbj, ubj, &
1207 & fscl, fmin, fmax, &
1209 & extract(ng) % rmask, &
1212 & extract(ng) % latr, &
1215 & extract(ng) % latr)
1217 IF (
founderror(status, nf90_noerr, __line__, myfile))
THEN
1223 WRITE (
stdout,30)
'latitude of RHO-points lat_rho', &
1224 & ng, trim(ncname), fmin, fmax
1230 extract(ng)%LatMin(ng)=fmin
1231 extract(ng)%LatMax(ng)=fmax
1233 CALL mp_exchange2d_xtr (ng, tile, model, 1, &
1234 & lbi, ubi, lbj, ubj, &
1236 & .false., .false., &
1237 & extract(ng) % latr)
1246 status=nf_fread2d_xtr(ng, model, ncname,
grx(ng)%ncid, &
1248 & 0, gtype, vsize, &
1249 & lbi, ubi, lbj, ubj, &
1250 & fscl, fmin, fmax, &
1252 & extract(ng) % umask, &
1255 & extract(ng) % lonu, &
1258 & extract(ng) % lonu)
1260 IF (
founderror(status, nf90_noerr, __line__, myfile))
THEN
1266 WRITE (
stdout,30)
'longitude of U-points: lon_u', &
1267 & ng, trim(ncname), fmin, fmax
1274 CALL mp_exchange2d_xtr (ng, tile, model, 1, &
1275 & lbi, ubi, lbj, ubj, &
1277 & .false., .false., &
1278 & extract(ng) % lonu)
1287 status=nf_fread2d_xtr(ng, model, ncname,
grx(ng)%ncid, &
1289 & 0, gtype, vsize, &
1290 & lbi, ubi, lbj, ubj, &
1291 & fscl, fmin, fmax, &
1293 & extract(ng) % umask, &
1296 & extract(ng) % latu, &
1299 & extract(ng) % latu)
1301 IF (
founderror(status, nf90_noerr, __line__, myfile))
THEN
1307 WRITE (
stdout,30)
'latitude of U-points: lat_u', &
1308 & ng, trim(ncname), fmin, fmax
1315 CALL mp_exchange2d_xtr (ng, tile, model, 1, &
1316 & lbi, ubi, lbj, ubj, &
1318 & .false., .false., &
1319 & extract(ng) % latu)
1328 status=nf_fread2d_xtr(ng, model, ncname,
grx(ng)%ncid, &
1330 & 0, gtype, vsize, &
1331 & lbi, ubi, lbj, ubj, &
1332 & fscl, fmin, fmax, &
1334 & extract(ng) % vmask, &
1337 & extract(ng) % lonv, &
1340 & extract(ng) % lonv)
1342 IF (
founderror(status, nf90_noerr, __line__, myfile))
THEN
1348 WRITE (
stdout,30)
'longitude of V-points: lon_v', &
1349 & ng, trim(ncname), fmin, fmax
1356 CALL mp_exchange2d_xtr (ng, tile, model, 1, &
1357 & lbi, ubi, lbj, ubj, &
1359 & .false., .false., &
1360 & extract(ng) % lonv)
1369 status=nf_fread2d_xtr(ng, model, ncname,
grx(ng)%ncid, &
1371 & 0, gtype, vsize, &
1372 & lbi, ubi, lbj, ubj, &
1373 & fscl, fmin, fmax, &
1375 & extract(ng) % vmask, &
1378 & extract(ng) % latv, &
1381 & extract(ng) % latv)
1383 IF (
founderror(status, nf90_noerr, __line__, myfile))
THEN
1389 WRITE (
stdout,30)
'latitude of V-points: lat_v', &
1390 & ng, trim(ncname), fmin, fmax
1397 CALL mp_exchange2d_xtr (ng, tile, model, 1, &
1398 & lbi, ubi, lbj, ubj, &
1400 & .false., .false., &
1401 & extract(ng) % latv)
1409 status=nf_fread2d_xtr(ng, model, ncname,
grx(ng)%ncid, &
1411 & 0, gtype, vsize, &
1412 & lbi, ubi, lbj, ubj, &
1413 & fscl, fmin, fmax, &
1415 & extract(ng) % rmask, &
1418 & extract(ng) % angler, &
1421 & extract(ng) % angler)
1423 IF (
founderror(status, nf90_noerr, __line__, myfile))
THEN
1429 WRITE (
stdout,30)
'angle between XI-axis and EAST: '// &
1431 & ng, trim(ncname), fmin, fmax
1437 IF (ewperiodic(ng).or.nsperiodic(ng))
THEN
1438 CALL exchange_r2d_xtr_tile (ng, tile, &
1439 & lbi, ubi, lbj, ubj, &
1440 & extract(ng) % angler)
1443 CALL mp_exchange2d_xtr (ng, tile, model, 1, &
1444 & lbi, ubi, lbj, ubj, &
1446 & ewperiodic(ng), nsperiodic(ng), &
1447 & extract(ng) % angler)
1455 IF (
founderror(exit_flag, noerror, __line__, myfile))
RETURN
1457 10
FORMAT (/,
' GET_EXTRACT_NF90 - unable to open grid NetCDF file:', &
1459 20
FORMAT (/,
' GET_EXTRACT_NF90 - unable to find grid variable: ', &
1460 & a,/,20x,
'in grid NetCDF file: ',a)
1461 30
FORMAT (2x,
'GET_EXTRACT_NF90 - ',a,/,22x, &
1462 &
'(Grid = ',i2.2,
', File: ',a,
')',/,22x, &
1463 &
'(Min = ', 1p,e15.8,0p,
' Max = ',1p,e15.8,0p,
')')
1464 40
FORMAT (/,
' GET_EXTRACT_NF90 - error while reading variable: ', &
1465 & a,/,20x,
'in grid NetCDF file: ',a)
1466 50
FORMAT (/,2x,
'GET_EXTRACT_NF90 - Reading adjoint sensitivity', &
1467 &
' scope arrays from file:',/22x,a,/)
1469 60
FORMAT (22x,
'(CheckSum = ',i0,
')')
1473 END SUBROUTINE get_extract_nf90
1475# if defined PIO_LIB && defined DISTRIBUTE
1478 SUBROUTINE get_extract_pio (ng, tile, model, &
1479 & LBi, UBi, LBj, UBj)
1484 integer,
intent(in) :: ng, tile, model
1485 integer,
intent(in) :: LBi, UBi, LBj, UBj
1489 integer :: cr, i, status, vindex
1492 integer(i8b) :: Fhash
1495 real(dp),
parameter :: Fscl = 1.0_dp
1497 real(r8) :: Fmax, Fmin
1499 character (len=256) :: ncname
1501 character (len=*),
parameter :: MyFile = &
1502 & __FILE__//
", get_extract_pio"
1504 TYPE (IO_desc_t),
pointer :: ioDesc
1505 TYPE (My_VarDesc) :: pioVar
1514 IF (
founderror(exit_flag, noerror, __line__, myfile))
RETURN
1519 IF (
grx(ng)%pioFile%fh.eq.-1)
THEN
1521 IF (
founderror(exit_flag, noerror, __line__, myfile))
THEN
1522 WRITE (
stdout,10) trim(ncname)
1530 & piofile =
grx(ng)%pioFile)
1531 IF (
founderror(exit_flag, noerror, __line__, myfile))
RETURN
1538 IF (
master)
WRITE (
stdout,20)
'spherical', trim(ncname)
1557# if (defined CURVGRID && defined UV_ADV)
1604 &
'spherical', spherical, &
1605 & piofile =
grx(ng)%pioFile)
1606 IF (
founderror(exit_flag, noerror, __line__, myfile))
RETURN
1626 SELECT CASE (trim(adjustl(
var_name(i))))
1633 IF (kind(extract(ng)%h).eq.8)
THEN
1634 piovar%dkind=pio_double
1637 piovar%dkind=pio_real
1641 status=nf_fread2d_xtr(ng, model, ncname,
grx(ng)%pioFile, &
1643 & 0, iodesc, vsize, &
1644 & lbi, ubi, lbj, ubj, &
1645 & fscl, fmin, fmax, &
1647 & extract(ng) % rmask, &
1650 & extract(ng) % h, &
1655 IF (
founderror(status, nf90_noerr, __line__, myfile))
THEN
1660# ifdef SINGLE_PRECISION
1661 extract(ng)%Hmin(ng)=real(fmin,dp)
1662 extract(ng)%Hmax(ng)=real(fmax,dp)
1664 extract(ng)%Hmin(ng)=fmin
1665 extract(ng)%Hmax(ng)=fmax
1668 WRITE (
stdout,30)
'bathymetry at RHO-points: h', &
1669 & ng, trim(ncname), fmin, fmax
1675 IF (ewperiodic(ng).or.nsperiodic(ng))
THEN
1676 CALL exchange_r2d_xtr_tile (ng, tile, &
1677 & lbi, ubi, lbj, ubj, &
1681 CALL mp_exchange2d_xtr (ng, tile, model, 1, &
1682 & lbi, ubi, lbj, ubj, &
1684 & ewperiodic(ng), nsperiodic(ng), &
1694 IF (kind(extract(ng)%rmask).eq.8)
THEN
1695 piovar%dkind=pio_double
1698 piovar%dkind=pio_real
1702 status=nf_fread2d_xtr(ng, model, ncname,
grx(ng)%pioFile, &
1704 & 0, iodesc, vsize, &
1705 & lbi, ubi, lbj, ubj, &
1706 & fscl, fmin, fmax, &
1707 & extract(ng) % rmask, &
1709 & extract(ng) % rmask, &
1712 & extract(ng) % rmask)
1714 IF (
founderror(status, nf90_noerr, __line__, myfile))
THEN
1720 WRITE (
stdout,30)
'mask on RHO-points: mask_rho', &
1721 & ng, trim(ncname), fmin, fmax
1727 IF (ewperiodic(ng).or.nsperiodic(ng))
THEN
1728 CALL exchange_r2d_xtr_tile (ng, tile, &
1729 & lbi, ubi, lbj, ubj, &
1730 & extract(ng) % rmask)
1733 CALL mp_exchange2d_xtr (ng, tile, model, 1, &
1734 & lbi, ubi, lbj, ubj, &
1736 & ewperiodic(ng), nsperiodic(ng), &
1737 & extract(ng) % rmask)
1745 IF (kind(extract(ng)%umask).eq.8)
THEN
1746 piovar%dkind=pio_double
1749 piovar%dkind=pio_real
1753 status=nf_fread2d_xtr(ng, model, ncname,
grx(ng)%pioFile, &
1755 & 0, iodesc, vsize, &
1756 & lbi, ubi, lbj, ubj, &
1757 & fscl, fmin, fmax, &
1758 & extract(ng) % umask, &
1760 & extract(ng) % umask, &
1763 & extract(ng) % umask)
1765 IF (
founderror(status, nf90_noerr, __line__, myfile))
THEN
1771 WRITE (
stdout,30)
'mask on U-points: mask_u', &
1772 & ng, trim(ncname), fmin, fmax
1778 IF (ewperiodic(ng).or.nsperiodic(ng))
THEN
1779 CALL exchange_u2d_xtr_tile (ng, tile, &
1780 & lbi, ubi, lbj, ubj, &
1781 & extract(ng) % umask)
1784 CALL mp_exchange2d_xtr (ng, tile, model, 1, &
1785 & lbi, ubi, lbj, ubj, &
1787 & ewperiodic(ng), nsperiodic(ng), &
1788 & extract(ng) % umask)
1796 IF (kind(extract(ng)%vmask).eq.8)
THEN
1797 piovar%dkind=pio_double
1800 piovar%dkind=pio_real
1804 status=nf_fread2d_xtr(ng, model, ncname,
grx(ng)%pioFile, &
1806 & 0, iodesc, vsize, &
1807 & lbi, ubi, lbj, ubj, &
1808 & fscl, fmin, fmax, &
1809 & extract(ng) % vmask, &
1811 & extract(ng) % vmask, &
1814 & extract(ng) % vmask)
1816 IF (
founderror(status, nf90_noerr, __line__, myfile))
THEN
1822 WRITE (
stdout,30)
'mask on V-points: mask_v', &
1823 & ng, trim(ncname), fmin, fmax
1829 IF (ewperiodic(ng).or.nsperiodic(ng))
THEN
1830 CALL exchange_v2d_xtr_tile (ng, tile, &
1831 & lbi, ubi, lbj, ubj, &
1832 & extract(ng) % vmask)
1835 CALL mp_exchange2d_xtr (ng, tile, model, 1, &
1836 & lbi, ubi, lbj, ubj, &
1838 & ewperiodic(ng), nsperiodic(ng), &
1839 & extract(ng) % vmask)
1847 IF (kind(extract(ng)%pmask).eq.8)
THEN
1848 piovar%dkind=pio_double
1851 piovar%dkind=pio_real
1855 status=nf_fread2d_xtr(ng, model, ncname,
grx(ng)%pioFile, &
1857 & 0, iodesc, vsize, &
1858 & lbi, ubi, lbj, ubj, &
1859 & fscl, fmin, fmax, &
1860 & extract(ng) % pmask, &
1862 & extract(ng) % pmask, &
1865 & extract(ng) % pmask)
1867 IF (
founderror(status, nf90_noerr, __line__, myfile))
THEN
1873 WRITE (
stdout,30)
'mask on PSI-points: mask_psi', &
1874 & ng, trim(ncname), fmin, fmax
1880 IF (ewperiodic(ng).or.nsperiodic(ng))
THEN
1881 CALL exchange_p2d_xtr_tile (ng, tile, &
1882 & lbi, ubi, lbj, ubj, &
1883 & extract(ng) % pmask)
1886 CALL mp_exchange2d_xtr (ng, tile, model, 1, &
1887 & lbi, ubi, lbj, ubj, &
1889 & ewperiodic(ng), nsperiodic(ng), &
1890 & extract(ng) % pmask)
1899 IF (kind(extract(ng)%pn).eq.8)
THEN
1900 piovar%dkind=pio_double
1903 piovar%dkind=pio_real
1907 status=nf_fread2d_xtr(ng, model, ncname,
grx(ng)%pioFile, &
1909 & 0, iodesc, vsize, &
1910 & lbi, ubi, lbj, ubj, &
1911 & fscl, fmin, fmax, &
1913 & extract(ng) % rmask, &
1916 & extract(ng) % f, &
1921 IF (
founderror(status, nf90_noerr, __line__, myfile))
THEN
1927 WRITE (
stdout,30)
'Coriolis parameter at RHO-points: f' &
1928 & ng, trim(ncname), fmin, fmax
1934 IF (ewperiodic(ng).or.nsperiodic(ng))
THEN
1935 CALL exchange_r2d_xtr_tile (ng, tile, &
1936 & lbi, ubi, lbj, ubj, &
1940 CALL mp_exchange2d_xtr (ng, tile, model, 1, &
1941 & lbi, ubi, lbj, ubj, &
1943 & ewperiodic(ng), nsperiodic(ng), &
1953 IF (kind(extract(ng)%pn).eq.8)
THEN
1954 piovar%dkind=pio_double
1957 piovar%dkind=pio_real
1961 status=nf_fread2d_xtr(ng, model, ncname,
grx(ng)%pioFile, &
1963 & 0, iodesc, vsize, &
1964 & lbi, ubi, lbj, ubj, &
1965 & fscl, fmin, fmax, &
1967 & extract(ng) % rmask, &
1970 & extract(ng) % pm, &
1975 IF (
founderror(status, nf90_noerr, __line__, myfile))
THEN
1981 WRITE (
stdout,30)
'reciprocal XI-grid spacing: pm', &
1982 & ng, trim(ncname), fmin, fmax
1988 IF (ewperiodic(ng).or.nsperiodic(ng))
THEN
1989 CALL exchange_r2d_xtr_tile (ng, tile, &
1990 & lbi, ubi, lbj, ubj, &
1994 CALL mp_exchange2d_xtr (ng, tile, model, 1, &
1995 & lbi, ubi, lbj, ubj, &
1997 & ewperiodic(ng), nsperiodic(ng), &
2007 IF (kind(extract(ng)%pn).eq.8)
THEN
2008 piovar%dkind=pio_double
2011 piovar%dkind=pio_real
2015 status=nf_fread2d_xtr(ng, model, ncname,
grx(ng)%pioFile, &
2017 & 0, iodesc, vsize, &
2018 & lbi, ubi, lbj, ubj, &
2019 & fscl, fmin, fmax, &
2021 & extract(ng) % rmask, &
2024 & extract(ng) % pn, &
2029 IF (
founderror(status, nf90_noerr, __line__, myfile))
THEN
2035 WRITE (
stdout,30)
'reciprocal ETA-grid spacing: pn', &
2036 & ng, trim(ncname), fmin, fmax
2042 IF (ewperiodic(ng).or.nsperiodic(ng))
THEN
2043 CALL exchange_r2d_xtr_tile (ng, tile, &
2044 & lbi, ubi, lbj, ubj, &
2048 CALL mp_exchange2d_xtr (ng, tile, model, 1, &
2049 & lbi, ubi, lbj, ubj, &
2051 & ewperiodic(ng), nsperiodic(ng), &
2054# if (defined CURVGRID && defined UV_ADV)
2061 IF (kind(extract(ng)%dmde).eq.8)
THEN
2062 piovar%dkind=pio_double
2065 piovar%dkind=pio_real
2069 status=nf_fread2d_xtr(ng, model, ncname,
grx(ng)%pioFile, &
2071 & 0, iodesc, vsize, &
2072 & lbi, ubi, lbj, ubj, &
2073 & fscl, fmin, fmax, &
2075 & extract(ng) % rmask, &
2078 & extract(ng) % dmde, &
2081 & extract(ng) % dmde)
2083 IF (
founderror(status, nf90_noerr, __line__, myfile))
THEN
2089 WRITE (
stdout,30)
'ETA-derivative of inverse metric '// &
2090 &
'factor pm: dmde', &
2091 & ng, trim(ncname), fmin, fmax
2097 IF (ewperiodic(ng).or.nsperiodic(ng))
THEN
2098 CALL exchange_r2d_xtr_tile (ng, tile, &
2099 & lbi, ubi, lbj, ubj, &
2100 & extract(ng) % dmde)
2103 CALL mp_exchange2d_xtr (ng, tile, model, 1, &
2104 & lbi, ubi, lbj, ubj, &
2106 & ewperiodic(ng), nsperiodic(ng), &
2107 & extract(ng) % dmde)
2115 IF (kind(extract(ng)%dndx).eq.8)
THEN
2116 piovar%dkind=pio_double
2119 piovar%dkind=pio_real
2123 status=nf_fread2d_xtr(ng, model, ncname,
grx(ng)%pioFile, &
2125 & 0, iodesc, vsize, &
2126 & lbi, ubi, lbj, ubj, &
2127 & fscl, fmin, fmax, &
2129 & extract(ng) % rmask, &
2132 & extract(ng) % dndx, &
2135 & extract(ng) % dndx)
2137 IF (
founderror(status, nf90_noerr, __line__, myfile))
THEN
2143 WRITE (
stdout,30)
'XI-derivative of inverse metric '// &
2144 &
'factor pn: dndx', &
2145 & ng, trim(ncname), fmin, fmax
2151 IF (ewperiodic(ng).or.nsperiodic(ng))
THEN
2152 CALL exchange_r2d_xtr_tile (ng, tile, &
2153 & lbi, ubi, lbj, ubj, &
2154 & extract(ng) % dndx)
2157 CALL mp_exchange2d_xtr (ng, tile, model, 1, &
2158 & lbi, ubi, lbj, ubj, &
2160 & ewperiodic(ng), nsperiodic(ng), &
2161 & extract(ng) % dndx)
2170 IF (kind(extract(ng)%xp).eq.8)
THEN
2171 piovar%dkind=pio_double
2174 piovar%dkind=pio_real
2178 status=nf_fread2d_xtr(ng, model, ncname,
grx(ng)%pioFile, &
2180 & 0, iodesc, vsize, &
2181 & lbi, ubi, lbj, ubj, &
2182 & fscl, fmin, fmax, &
2184 & extract(ng) % pmask, &
2187 & extract(ng) % xp, &
2192 IF (
founderror(status, nf90_noerr, __line__, myfile))
THEN
2198 WRITE (
stdout,30)
'x-location of PSI-points: x_psi', &
2199 & ng, trim(ncname), fmin, fmax
2206 CALL mp_exchange2d_xtr (ng, tile, model, 1, &
2207 & lbi, ubi, lbj, ubj, &
2209 & .false., .false., &
2218 IF (kind(extract(ng)%yp).eq.8)
THEN
2219 piovar%dkind=pio_double
2222 piovar%dkind=pio_real
2226 status=nf_fread2d_xtr(ng, model, ncname,
grx(ng)%pioFile, &
2228 & 0, iodesc, vsize, &
2229 & lbi, ubi, lbj, ubj, &
2230 & fscl, fmin, fmax, &
2232 & extract(ng) % pmask, &
2235 & extract(ng) % yp, &
2240 IF (
founderror(status, nf90_noerr, __line__, myfile))
THEN
2246 WRITE (
stdout,30)
'y-location of PSI-points: y-psi', &
2247 & ng, trim(ncname), fmin, fmax
2254 CALL mp_exchange2d_xtr (ng, tile, model, 1, &
2255 & lbi, ubi, lbj, ubj, &
2257 & .false., .false., &
2266 IF (kind(extract(ng)%xr).eq.8)
THEN
2267 piovar%dkind=pio_double
2270 piovar%dkind=pio_real
2274 status=nf_fread2d_xtr(ng, model, ncname,
grx(ng)%pioFile, &
2276 & 0, iodesc, vsize, &
2277 & lbi, ubi, lbj, ubj, &
2278 & fscl, fmin, fmax, &
2280 & extract(ng) % rmask, &
2283 & extract(ng) % xr, &
2288 IF (
founderror(status, nf90_noerr, __line__, myfile))
THEN
2294 WRITE (
stdout,30)
'x-location of RHO-points: x-rho', &
2295 & ng, trim(ncname), fmin, fmax
2301 IF (.not.spherical)
THEN
2302 extract(ng)%LonMin(ng)=fmin
2303 extract(ng)%LonMax(ng)=fmax
2306 CALL mp_exchange2d_xtr (ng, tile, model, 1, &
2307 & lbi, ubi, lbj, ubj, &
2309 & .false., .false., &
2318 IF (kind(extract(ng)%yr).eq.8)
THEN
2319 piovar%dkind=pio_double
2322 piovar%dkind=pio_real
2326 status=nf_fread2d_xtr(ng, model, ncname,
grx(ng)%pioFile, &
2328 & 0, iodesc, vsize, &
2329 & lbi, ubi, lbj, ubj, &
2330 & fscl, fmin, fmax, &
2332 & extract(ng) % rmask, &
2335 & extract(ng) % yr, &
2340 IF (
founderror(status, nf90_noerr, __line__, myfile))
THEN
2346 WRITE (
stdout,30)
'y-location of RHO-points: y_rho', &
2347 & ng, trim(ncname), fmin, fmax
2353 IF (.not.spherical)
THEN
2354 extract(ng)%LatMin(ng)=fmin
2355 extract(ng)%LatMax(ng)=fmax
2358 CALL mp_exchange2d_xtr (ng, tile, model, 1, &
2359 & lbi, ubi, lbj, ubj, &
2361 & .false., .false., &
2370 IF (kind(extract(ng)%xu).eq.8)
THEN
2371 piovar%dkind=pio_double
2374 piovar%dkind=pio_real
2378 status=nf_fread2d_xtr(ng, model, ncname,
grx(ng)%pioFile, &
2380 & 0, iodesc, vsize, &
2381 & lbi, ubi, lbj, ubj, &
2382 & fscl, fmin, fmax, &
2384 & extract(ng) % umask, &
2387 & extract(ng) % xu, &
2392 IF (
founderror(status, nf90_noerr, __line__, myfile))
THEN
2398 WRITE (
stdout,30)
'x-location of U-points: x_u', &
2399 & ng, trim(ncname), fmin, fmax
2406 CALL mp_exchange2d_xtr (ng, tile, model, 1, &
2407 & lbi, ubi, lbj, ubj, &
2409 & .false., .false., &
2418 IF (kind(extract(ng)%yu).eq.8)
THEN
2419 piovar%dkind=pio_double
2422 piovar%dkind=pio_real
2426 status=nf_fread2d_xtr(ng, model, ncname,
grx(ng)%pioFile, &
2428 & 0, iodesc, vsize, &
2429 & lbi, ubi, lbj, ubj, &
2430 & fscl, fmin, fmax, &
2432 & extract(ng) % umask, &
2435 & extract(ng) % yu, &
2440 IF (
founderror(status, nf90_noerr, __line__, myfile))
THEN
2446 WRITE (
stdout,30)
'y-location of U-points: y_u', &
2447 & ng, trim(ncname), fmin, fmax
2454 CALL mp_exchange2d_xtr (ng, tile, model, 1, &
2455 & lbi, ubi, lbj, ubj, &
2457 & .false., .false., &
2466 IF (kind(extract(ng)%xv).eq.8)
THEN
2467 piovar%dkind=pio_double
2470 piovar%dkind=pio_real
2474 status=nf_fread2d_xtr(ng, model, ncname,
grx(ng)%pioFile, &
2476 & 0, iodesc, vsize, &
2477 & lbi, ubi, lbj, ubj, &
2478 & fscl, fmin, fmax, &
2480 & extract(ng) % vmask, &
2483 & extract(ng) % xv, &
2488 IF (
founderror(status, nf90_noerr, __line__, myfile))
THEN
2494 WRITE (
stdout,30)
'x-location of V-points: x_v', &
2495 & ng, trim(ncname), fmin, fmax
2502 CALL mp_exchange2d_xtr (ng, tile, model, 1, &
2503 & lbi, ubi, lbj, ubj, &
2505 & .false., .false., &
2514 IF (kind(extract(ng)%yv).eq.8)
THEN
2515 piovar%dkind=pio_double
2518 piovar%dkind=pio_real
2522 status=nf_fread2d_xtr(ng, model, ncname,
grx(ng)%pioFile, &
2524 & 0, iodesc, vsize, &
2525 & lbi, ubi, lbj, ubj, &
2526 & fscl, fmin, fmax, &
2528 & extract(ng) % vmask, &
2531 & extract(ng) % yv, &
2536 IF (
founderror(status, nf90_noerr, __line__, myfile))
THEN
2542 WRITE (
stdout,30)
'y-location of V-points: y_v', &
2543 & ng, trim(ncname), fmin, fmax
2550 CALL mp_exchange2d_xtr (ng, tile, model, 1, &
2551 & lbi, ubi, lbj, ubj, &
2553 & .false., .false., &
2563 IF (kind(extract(ng)%lonp).eq.8)
THEN
2564 piovar%dkind=pio_double
2567 piovar%dkind=pio_real
2571 status=nf_fread2d_xtr(ng, model, ncname,
grx(ng)%pioFile, &
2573 & 0, iodesc, vsize, &
2574 & lbi, ubi, lbj, ubj, &
2575 & fscl, fmin, fmax, &
2577 & extract(ng) % pmask, &
2580 & extract(ng) % lonp, &
2583 & extract(ng) % lonp)
2585 IF (
founderror(status, nf90_noerr, __line__, myfile))
THEN
2591 WRITE (
stdout,30)
'longitude of PSI-points: lon_psi', &
2592 & ng, trim(ncname), fmin, fmax
2599 CALL mp_exchange2d_xtr (ng, tile, model, 1, &
2600 & lbi, ubi, lbj, ubj, &
2602 & .false., .false., &
2603 & extract(ng) % lonp)
2613 IF (kind(extract(ng)%latp).eq.8)
THEN
2614 piovar%dkind=pio_double
2617 piovar%dkind=pio_real
2621 status=nf_fread2d_xtr(ng, model, ncname,
grx(ng)%pioFile, &
2623 & 0, iodesc, vsize, &
2624 & lbi, ubi, lbj, ubj, &
2625 & fscl, fmin, fmax, &
2627 & extract(ng) % pmask, &
2630 & extract(ng) % latp, &
2633 & extract(ng) % latp)
2635 IF (
founderror(status, nf90_noerr, __line__, myfile))
THEN
2641 WRITE (
stdout,30)
'latitude of PSI-points lat_psi', &
2642 & ng, trim(ncname), fmin, fmax
2649 CALL mp_exchange2d_xtr (ng, tile, model, 1, &
2650 & lbi, ubi, lbj, ubj, &
2652 & .false., .false., &
2653 & extract(ng) % latp)
2663 IF (kind(extract(ng)%lonr).eq.8)
THEN
2664 piovar%dkind=pio_double
2667 piovar%dkind=pio_real
2671 status=nf_fread2d_xtr(ng, model, ncname,
grx(ng)%pioFile, &
2673 & 0, iodesc, vsize, &
2674 & lbi, ubi, lbj, ubj, &
2675 & fscl, fmin, fmax, &
2677 & extract(ng) % rmask, &
2680 & extract(ng) % lonr, &
2683 & extract(ng) % lonr)
2685 IF (
founderror(status, nf90_noerr, __line__, myfile))
THEN
2691 WRITE (
stdout,30)
'longitude of RHO-points: lon_rho', &
2692 & ng, trim(ncname), fmin, fmax
2698 extract(ng)%LonMin(ng)=fmin
2699 extract(ng)%LonMax(ng)=fmax
2701 CALL mp_exchange2d_xtr (ng, tile, model, 1, &
2702 & lbi, ubi, lbj, ubj, &
2704 & .false., .false., &
2705 & extract(ng) % lonr)
2715 IF (kind(extract(ng)%latr).eq.8)
THEN
2716 piovar%dkind=pio_double
2719 piovar%dkind=pio_real
2723 status=nf_fread2d_xtr(ng, model, ncname,
grx(ng)%pioFile, &
2725 & 0, iodesc, vsize, &
2726 & lbi, ubi, lbj, ubj, &
2727 & fscl, fmin, fmax, &
2729 & extract(ng) % rmask, &
2732 & extract(ng) % latr, &
2735 & extract(ng) % latr)
2737 IF (
founderror(status, nf90_noerr, __line__, myfile))
THEN
2743 WRITE (
stdout,30)
'latitude of RHO-points lat_rho', &
2744 & ng, trim(ncname), fmin, fmax
2750 extract(ng)%LatMin(ng)=fmin
2751 extract(ng)%LatMax(ng)=fmax
2753 CALL mp_exchange2d_xtr (ng, tile, model, 1, &
2754 & lbi, ubi, lbj, ubj, &
2756 & .false., .false., &
2757 & extract(ng) % latr)
2767 IF (kind(extract(ng)%lonu).eq.8)
THEN
2768 piovar%dkind=pio_double
2771 piovar%dkind=pio_real
2775 status=nf_fread2d_xtr(ng, model, ncname,
grx(ng)%pioFile, &
2777 & 0, iodesc, vsize, &
2778 & lbi, ubi, lbj, ubj, &
2779 & fscl, fmin, fmax, &
2781 & extract(ng) % umask, &
2784 & extract(ng) % lonu, &
2787 & extract(ng) % lonu)
2789 IF (
founderror(status, nf90_noerr, __line__, myfile))
THEN
2795 WRITE (
stdout,30)
'longitude of U-points: lon_u', &
2796 & ng, trim(ncname), fmin, fmax
2803 CALL mp_exchange2d_xtr (ng, tile, model, 1, &
2804 & lbi, ubi, lbj, ubj, &
2806 & .false., .false., &
2807 & extract(ng) % lonu)
2817 IF (kind(extract(ng)%latu).eq.8)
THEN
2818 piovar%dkind=pio_double
2821 piovar%dkind=pio_real
2825 status=nf_fread2d_xtr(ng, model, ncname,
grx(ng)%pioFile, &
2827 & 0, iodesc, vsize, &
2828 & lbi, ubi, lbj, ubj, &
2829 & fscl, fmin, fmax, &
2831 & extract(ng) % umask, &
2834 & extract(ng) % latu, &
2837 & extract(ng) % latu)
2839 IF (
founderror(status, nf90_noerr, __line__, myfile))
THEN
2845 WRITE (
stdout,30)
'latitude of U-points: lat_u', &
2846 & ng, trim(ncname), fmin, fmax
2853 CALL mp_exchange2d_xtr (ng, tile, model, 1, &
2854 & lbi, ubi, lbj, ubj, &
2856 & .false., .false., &
2857 & extract(ng) % latu)
2867 IF (kind(extract(ng)%lonv).eq.8)
THEN
2868 piovar%dkind=pio_double
2871 piovar%dkind=pio_real
2875 status=nf_fread2d_xtr(ng, model, ncname,
grx(ng)%pioFile, &
2877 & 0, iodesc, vsize, &
2878 & lbi, ubi, lbj, ubj, &
2879 & fscl, fmin, fmax, &
2881 & extract(ng) % vmask, &
2884 & extract(ng) % lonv, &
2887 & extract(ng) % lonv)
2889 IF (
founderror(status, nf90_noerr, __line__, myfile))
THEN
2895 WRITE (
stdout,30)
'longitude of V-points: lon_v', &
2896 & ng, trim(ncname), fmin, fmax
2903 CALL mp_exchange2d_xtr (ng, tile, model, 1, &
2904 & lbi, ubi, lbj, ubj, &
2906 & .false., .false., &
2907 & extract(ng) % lonv)
2917 IF (kind(extract(ng)%latv).eq.8)
THEN
2918 piovar%dkind=pio_double
2921 piovar%dkind=pio_real
2925 status=nf_fread2d_xtr(ng, model, ncname,
grx(ng)%pioFile, &
2927 & 0, iodesc, vsize, &
2928 & lbi, ubi, lbj, ubj, &
2929 & fscl, fmin, fmax, &
2931 & extract(ng) % vmask, &
2934 & extract(ng) % latv, &
2937 & extract(ng) % latv)
2939 IF (
founderror(status, nf90_noerr, __line__, myfile))
THEN
2945 WRITE (
stdout,30)
'latitude of V-points: lat_v', &
2946 & ng, trim(ncname), fmin, fmax
2953 CALL mp_exchange2d_xtr (ng, tile, model, 1, &
2954 & lbi, ubi, lbj, ubj, &
2956 & .false., .false., &
2957 & extract(ng) % latv)
2966 IF (kind(extract(ng)%angler).eq.8)
THEN
2967 piovar%dkind=pio_double
2970 piovar%dkind=pio_real
2974 status=nf_fread2d_xtr(ng, model, ncname,
grx(ng)%pioFile, &
2976 & 0, iodesc, vsize, &
2977 & lbi, ubi, lbj, ubj, &
2978 & fscl, fmin, fmax, &
2980 & extract(ng) % rmask, &
2983 & extract(ng) % angler, &
2986 & extract(ng) % angler)
2988 IF (
founderror(status, nf90_noerr, __line__, myfile))
THEN
2994 WRITE (
stdout,30)
'angle between XI-axis and EAST: '// &
2996 & ng, trim(ncname), fmin, fmax
3002 IF (ewperiodic(ng).or.nsperiodic(ng))
THEN
3003 CALL exchange_r2d_xtr_tile (ng, tile, &
3004 & lbi, ubi, lbj, ubj, &
3005 & extract(ng) % angler)
3008 CALL mp_exchange2d_xtr (ng, tile, model, 1, &
3009 & lbi, ubi, lbj, ubj, &
3011 & ewperiodic(ng), nsperiodic(ng), &
3012 & extract(ng) % angler)
3020 IF (
founderror(exit_flag, noerror, __line__, myfile))
RETURN
3022 10
FORMAT (/,
' GET_EXTRACT_PIO - unable to open grid NetCDF file:', &
3024 20
FORMAT (/,
' GET_EXTRACT_PIO - unable to find grid variable: ',a, &
3025 & /,19x,
'in grid NetCDF file: ',a)
3026 30
FORMAT (2x,
'GET_EXTRACT_PIO - ',a,/,22x, &
3027 &
'(Grid = ',i2.2,
', File: ',a,
')',/,22x, &
3028 &
'(Min = ', 1p,e15.8,0p,
' Max = ',1p,e15.8,0p,
')')
3029 40
FORMAT (/,
' GET_EXTRACT_PIO - error while reading variable: ',a, &
3030 & /,12x,
'in grid NetCDF file: ',a)
3031 50
FORMAT (/,2x,
'GET_EXTRACT_PIO - Reading adjoint sensitivity', &
3032 &
' scope arrays from file:',/22x,a,/)
3034 60
FORMAT (22x,
'(CheckSum = ',i0,
')')
3038 END SUBROUTINE get_extract_pio
type(t_grid), dimension(:), allocatable grid
type(t_io), dimension(:), allocatable grx
character(len=256) sourcefile
integer, parameter io_nf90
integer, parameter io_pio
subroutine, public netcdf_close(ng, model, ncid, ncname, lupdate)
subroutine, public netcdf_open(ng, model, ncname, omode, ncid)
character(len=100), dimension(mvars) var_name
integer, dimension(mvars) var_id
subroutine, public netcdf_inq_var(ng, model, ncname, ncid, myvarname, searchvar, varid, nvardim, nvaratt)
integer, parameter u2dvar
integer, parameter p2dvar
integer, parameter r2dvar
integer, parameter v2dvar
type(io_desc_t), dimension(:), pointer iodesc_sp_u2dvar
type(var_desc_t), dimension(:), pointer var_desc
type(io_desc_t), dimension(:), pointer iodesc_sp_r2dvar
subroutine, public pio_netcdf_inq_var(ng, model, ncname, piofile, myvarname, searchvar, piovar, nvardim, nvaratt)
subroutine, public pio_netcdf_close(ng, model, piofile, ncname, lupdate)
subroutine, public pio_netcdf_open(ng, model, ncname, omode, piofile)
type(io_desc_t), dimension(:), pointer iodesc_dp_u2dvar
type(io_desc_t), dimension(:), pointer iodesc_dp_p2dvar
type(io_desc_t), dimension(:), pointer iodesc_dp_v2dvar
type(io_desc_t), dimension(:), pointer iodesc_dp_r2dvar
type(io_desc_t), dimension(:), pointer iodesc_sp_p2dvar
type(io_desc_t), dimension(:), pointer iodesc_sp_v2dvar
logical function, public find_string(a, asize, string, aindex)
logical function, public founderror(flag, noerr, line, routine)