Changeset 2408 for LMDZ5/branches/testing/libf/phylmd/regr_pr_o3_m.F90
- Timestamp:
- Dec 14, 2015, 11:43:09 AM (9 years ago)
- Location:
- LMDZ5/branches/testing
- Files:
-
- 2 edited
Legend:
- Unmodified
- Added
- Removed
-
LMDZ5/branches/testing
- Property svn:mergeinfo changed
/LMDZ5/trunk merged: 2293-2295,2297,2299-2302,2305-2313,2315,2317-2380,2382-2396
- Property svn:mergeinfo changed
-
LMDZ5/branches/testing/libf/phylmd/regr_pr_o3_m.F90
r1910 r2408 30 30 use regr1_step_av_m, only: regr1_step_av 31 31 use press_coefoz_m, only: press_in_edg 32 use control_mod, only: dayref 32 use time_phylmdz_mod, only: day_ref 33 use mod_grid_phy_lmdz, only: nbp_lon, nbp_lat, nbp_lev 33 34 34 35 REAL, intent(in):: p3d(:, :, :) ! pressure at layer interfaces, in Pa … … 43 44 ! Variables local to the procedure: 44 45 45 include "dimensions.h"46 47 46 integer ncid, varid, ncerr ! for NetCDF 48 47 integer i, j 49 48 50 real r_mob( jjm + 1, size(press_in_edg) - 1)51 ! (ozone mole fraction from Mobidic at day "day ref")49 real r_mob(nbp_lat, size(press_in_edg) - 1) 50 ! (ozone mole fraction from Mobidic at day "day_ref") 52 51 ! (r_mob(j, k) is at latitude "rlatu(j)", in pressure interval 53 52 ! "[press_in_edg(k), press_in_edg(k+1)]".) … … 56 55 57 56 print *, "Call sequence information: regr_pr_o3" 58 call assert(shape(o3_mob_regr) == (/ iim + 1, jjm + 1, llm/), &57 call assert(shape(o3_mob_regr) == (/nbp_lon + 1, nbp_lat, nbp_lev/), & 59 58 "regr_pr_o3 o3_mob_regr") 60 call assert(shape(p3d) == (/ iim + 1, jjm + 1, llm+ 1/), "regr_pr_o3 p3d")59 call assert(shape(p3d) == (/nbp_lon + 1, nbp_lat, nbp_lev + 1/), "regr_pr_o3 p3d") 61 60 62 61 call nf95_open("coefoz_LMDZ.nc", nf90_nowrite, ncid) … … 64 63 call nf95_inq_varid(ncid, "r_Mob", varid) 65 64 ! Get data at the right day from the input file: 66 ncerr = nf90_get_var(ncid, varid, r_mob, start=(/1, 1, day ref/))65 ncerr = nf90_get_var(ncid, varid, r_mob, start=(/1, 1, day_ref/)) 67 66 call handle_err("nf90_get_var r_Mob", ncerr) 68 67 ! Latitudes are in ascending order in the input file while 69 68 ! "rlatu" is in descending order so we need to invert order: 70 r_mob = r_mob( jjm+1:1:-1, :)69 r_mob = r_mob(nbp_lat:1:-1, :) 71 70 72 71 call nf95_close(ncid) … … 75 74 76 75 ! Poles: 77 do j = 1, jjm + 1, jjm78 o3_mob_regr(1, j, llm:1:-1) &79 = regr1_step_av(r_mob(j, :), press_in_edg, p3d(1, j, llm+1:1:-1))76 do j = 1, nbp_lat, nbp_lat-1 77 o3_mob_regr(1, j, nbp_lev:1:-1) & 78 = regr1_step_av(r_mob(j, :), press_in_edg, p3d(1, j, nbp_lev+1:1:-1)) 80 79 ! (invert order of indices because "p3d" is in descending order) 81 80 end do 82 81 83 82 ! Other latitudes: 84 do j = 2, jjm85 do i = 1, iim86 o3_mob_regr(i, j, llm:1:-1) &83 do j = 2, nbp_lat-1 84 do i = 1, nbp_lon 85 o3_mob_regr(i, j, nbp_lev:1:-1) & 87 86 = regr1_step_av(r_mob(j, :), press_in_edg, & 88 p3d(i, j, llm+1:1:-1))87 p3d(i, j, nbp_lev+1:1:-1)) 89 88 ! (invert order of indices because "p3d" is in descending order) 90 89 end do … … 92 91 93 92 ! Duplicate pole values on all longitudes: 94 o3_mob_regr(2:, 1, :) = spread(o3_mob_regr(1, 1, :), dim=1, ncopies= iim)95 o3_mob_regr(2:, jjm + 1, :) &96 = spread(o3_mob_regr(1, jjm + 1, :), dim=1, ncopies=iim)93 o3_mob_regr(2:, 1, :) = spread(o3_mob_regr(1, 1, :), dim=1, ncopies=nbp_lon) 94 o3_mob_regr(2:, nbp_lat, :) & 95 = spread(o3_mob_regr(1, nbp_lat, :), dim=1, ncopies=nbp_lon) 97 96 98 97 ! Duplicate first longitude to last longitude: 99 o3_mob_regr( iim + 1, 2:jjm, :) = o3_mob_regr(1, 2:jjm, :)98 o3_mob_regr(nbp_lon + 1, 2:nbp_lat-1, :) = o3_mob_regr(1, 2:nbp_lat-1, :) 100 99 101 100 end subroutine regr_pr_o3
Note: See TracChangeset
for help on using the changeset viewer.