1 | ! |
---|
2 | ! $Id: stratH2O_methox.F90 3677 2020-05-06 15:18:32Z oboucher $ |
---|
3 | ! |
---|
4 | SUBROUTINE stratH2O_methox(debutphy,paprs,dq_ch4mmr) |
---|
5 | ! |
---|
6 | ! output: CH4VMR in mmr (mass mixing ratio/sec: kg H2O/kg air) |
---|
7 | |
---|
8 | USE netcdf95, ONLY: nf95_close, nf95_gw_var, nf95_inq_dimid, & |
---|
9 | nf95_inq_varid, nf95_inquire_dimension, nf95_open |
---|
10 | USE netcdf, ONLY: nf90_get_var, nf90_noerr, nf90_nowrite |
---|
11 | |
---|
12 | USE mod_grid_phy_lmdz |
---|
13 | USE mod_phys_lmdz_mpi_data, ONLY : is_mpi_root |
---|
14 | USE mod_phys_lmdz_omp_data, ONLY : is_omp_root |
---|
15 | |
---|
16 | USE mod_phys_lmdz_para |
---|
17 | USE dimphy |
---|
18 | USE phys_cal_mod, ONLY : mth_cur |
---|
19 | USE infotrac_phy |
---|
20 | USE aerophys |
---|
21 | USE YOMCST |
---|
22 | USE strataer_local_var_mod, ONLY : flag_newclim_file |
---|
23 | |
---|
24 | IMPLICIT NONE |
---|
25 | |
---|
26 | include "dimensions.h" |
---|
27 | |
---|
28 | ! Variable input |
---|
29 | REAL paprs(klon,klev+1) |
---|
30 | LOGICAL, INTENT(IN) :: debutphy ! le flag de l'initialisation de la physique |
---|
31 | ! Variable output |
---|
32 | ! tendance buffer pour appel de add_phys_tend |
---|
33 | REAL, INTENT(INOUT), DIMENSION(klon,klev) :: dq_ch4mmr |
---|
34 | |
---|
35 | ! Variables locales |
---|
36 | INTEGER n_lat ! number of latitudes in the input data |
---|
37 | INTEGER n_lon ! number of longitudes in the input data |
---|
38 | INTEGER, SAVE :: n_lev ! number of levels in the input data |
---|
39 | !$OMP THREADPRIVATE(n_lev) |
---|
40 | INTEGER n_mth ! number of months in the input data |
---|
41 | INTEGER, SAVE :: mth_pre |
---|
42 | !$OMP THREADPRIVATE(mth_pre) |
---|
43 | |
---|
44 | ! Champs reconstitues |
---|
45 | REAL paprs_glo(klon_glo,klev+1) |
---|
46 | |
---|
47 | REAL, ALLOCATABLE:: latitude(:) |
---|
48 | ! (of input data sorted in strictly ascending order) |
---|
49 | REAL, ALLOCATABLE:: longitude(:) |
---|
50 | ! (of input data sorted in strictly ascending order) |
---|
51 | REAL, ALLOCATABLE:: time(:) |
---|
52 | ! (of input data sorted in strictly ascending order) |
---|
53 | REAL, ALLOCATABLE:: lev(:) |
---|
54 | ! levels of input data |
---|
55 | |
---|
56 | !stratospheric H2O source from CH4 oxidation |
---|
57 | ! fixed climatos |
---|
58 | ! H2O production in VMR/sec) |
---|
59 | REAL, ALLOCATABLE :: CH4RVMR_in(:, :, :, :) |
---|
60 | REAL, ALLOCATABLE :: CH4RVMR_mth(:, :, :) |
---|
61 | REAL, ALLOCATABLE :: CH4RVMR_tmp(:, :) |
---|
62 | REAL CH4RVMR_glo(klon_glo,klev) |
---|
63 | ! |
---|
64 | INTEGER i, k, kk, j |
---|
65 | |
---|
66 | ! For NetCDF: |
---|
67 | INTEGER ncid_in ! IDs for input files |
---|
68 | INTEGER varid, ncerr |
---|
69 | |
---|
70 | INTEGER, PARAMETER :: lev_input=17 |
---|
71 | !--pressure at interfaces of input data (in Pa) |
---|
72 | REAL, DIMENSION(lev_input+1), PARAMETER :: & |
---|
73 | paprs_input=(/ & |
---|
74 | 1.00000002e+05, 6.06530673e+04, 3.67879449e+04, & |
---|
75 | 2.23130165e+04, 1.35335286e+04, 8.20850004e+03, & |
---|
76 | 4.97870695e+03, 3.01973841e+03, 1.83156393e+03, & |
---|
77 | 1.11089968e+03, 6.73794715e+02, 4.08677153e+02, & |
---|
78 | 2.47875223e+02, 1.50343923e+02, 9.11881985e+01, & |
---|
79 | 5.53084382e+01, 3.35462635e+01, 0.0 /) |
---|
80 | ! |
---|
81 | |
---|
82 | IF (debutphy .OR. mth_cur .NE. mth_pre) THEN |
---|
83 | |
---|
84 | !--preparation of global fields |
---|
85 | CALL gather(paprs, paprs_glo) |
---|
86 | |
---|
87 | IF (is_mpi_root.AND.is_omp_root) THEN |
---|
88 | |
---|
89 | print *,'In stratH2O_methox read file: ch4r_annual_lmdz.nc' |
---|
90 | |
---|
91 | !--reading strat. H2O source (CH4 oxidation) files |
---|
92 | CALL nf95_open("ch4r_annual_lmdz.nc", nf90_nowrite, ncid_in) |
---|
93 | |
---|
94 | CALL nf95_inq_varid(ncid_in, "LEV", varid) |
---|
95 | CALL nf95_gw_var(ncid_in, varid, lev) |
---|
96 | n_lev = size(lev) |
---|
97 | |
---|
98 | CALL nf95_inq_varid(ncid_in, "LAT", varid) |
---|
99 | CALL nf95_gw_var(ncid_in, varid, latitude) |
---|
100 | n_lat = size(latitude) |
---|
101 | |
---|
102 | CALL nf95_inq_varid(ncid_in, "LON", varid) |
---|
103 | CALL nf95_gw_var(ncid_in, varid, longitude) |
---|
104 | n_lon = size(longitude) |
---|
105 | |
---|
106 | CALL nf95_inq_varid(ncid_in, "TIME", varid) |
---|
107 | CALL nf95_gw_var(ncid_in, varid, time) |
---|
108 | n_mth = size(time) |
---|
109 | |
---|
110 | IF (.NOT.ALLOCATED(CH4RVMR_in)) ALLOCATE(CH4RVMR_in(n_lon, n_lat, n_lev, n_mth)) |
---|
111 | |
---|
112 | CALL nf95_inq_varid(ncid_in, "CH4RVMR", varid) |
---|
113 | ncerr = nf90_get_var(ncid_in, varid, CH4RVMR_in) |
---|
114 | print *,'code erreur CH4RVMR=', ncerr, varid |
---|
115 | |
---|
116 | CALL nf95_close(ncid_in) |
---|
117 | |
---|
118 | IF (.NOT.ALLOCATED(CH4RVMR_mth)) ALLOCATE(CH4RVMR_mth(n_lon, n_lat, n_lev)) |
---|
119 | IF (.NOT.ALLOCATED(CH4RVMR_tmp)) ALLOCATE(CH4RVMR_tmp(klon_glo, n_lev)) |
---|
120 | |
---|
121 | !---select the correct month, |
---|
122 | !---correct latitudinal order,convert input from volume mixing ratio to mass mixing ratio |
---|
123 | DO j=1,n_lat |
---|
124 | ! convert VMR/sec in mmr (mass mixing ratio/sec: kg H2O/kg air) |
---|
125 | ! x2 because CH4->2*H2O |
---|
126 | CH4RVMR_mth(:,j,:) = 2*CH4RVMR_in(:,n_lat+1-j,:,mth_cur)*mH2Omol/mAIRmol |
---|
127 | ENDDO |
---|
128 | |
---|
129 | !---reduce to a klon_glo grid but keep the levels |
---|
130 | CALL grid2dTo1d_glo(CH4RVMR_mth,CH4RVMR_tmp) |
---|
131 | |
---|
132 | !---regrid weighted climatologies |
---|
133 | DO i=1, klon_glo |
---|
134 | DO k=1, klev |
---|
135 | |
---|
136 | CH4RVMR_glo(i,k)=0.0 |
---|
137 | |
---|
138 | DO kk=1, n_lev |
---|
139 | ! |
---|
140 | CH4RVMR_glo(i,k)=CH4RVMR_glo(i,k)+ & |
---|
141 | MAX(0.0,MIN(paprs_glo(i,k),paprs_input(kk)) & |
---|
142 | -MAX(paprs_glo(i,k+1),paprs_input(kk+1))) & |
---|
143 | *CH4RVMR_tmp(i,kk)/(paprs_glo(i,k)-paprs_glo(i,k+1)) |
---|
144 | ENDDO ! kk loop |
---|
145 | ENDDO ! k loop |
---|
146 | ENDDO ! i loop |
---|
147 | ENDIF !--is_mpi_root and is_omp_root |
---|
148 | |
---|
149 | !--keep memory of previous month |
---|
150 | mth_pre=mth_cur |
---|
151 | |
---|
152 | !--scatter global fields around |
---|
153 | CALL scatter(CH4RVMR_glo, dq_ch4mmr) |
---|
154 | |
---|
155 | IF (is_mpi_root.AND.is_omp_root) THEN |
---|
156 | DEALLOCATE(CH4RVMR_in,CH4RVMR_mth,CH4RVMR_tmp) |
---|
157 | ENDIF !--is_mpi_root and is_omp_root |
---|
158 | ENDIF ! debutphy.OR.new month |
---|
159 | |
---|
160 | RETURN |
---|
161 | |
---|
162 | END SUBROUTINE stratH2O_methox |
---|