source: LMDZ6/trunk/libf/phylmd/StratAer/stratH2O_methox.f90 @ 5445

Last change on this file since 5445 was 5338, checked in by Laurent Fairhead, 5 weeks ago

Getting rid of dependance to the dynamics

File size: 5.4 KB
Line 
1!
2! $Id: stratH2O_methox.F90 3677 2020-05-06 15:18:32Z oboucher $
3!
4SUBROUTINE stratH2O_methox(debutphy,paprs,dq_ch4mmr)
5!
6! output: CH4VMR in MMR/s (mass mixing ratio/s or kg H2O/kg air/s)
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_mod_h
22  USE strataer_local_var_mod, ONLY : flag_newclim_file
23 
24IMPLICIT NONE
25
26
27 
28! Input variables
29  REAL paprs(klon,klev+1)
30  LOGICAL, INTENT(IN) :: debutphy   ! flag for first physiq step
31! Output variables
32! tendency buffer used in add_phys_tend subroutine (in physiq_mod)
33  REAL, INTENT(INOUT), DIMENSION(klon,klev)  :: dq_ch4mmr
34 
35! Local variables
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! Reconstitutes fields
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 from fixed climatos
57! (H2O production in VMR/sec)
58  REAL, ALLOCATABLE :: CH4RVMR_in(:, :, :, :)
59  REAL, ALLOCATABLE :: CH4RVMR_mth(:, :, :)
60  REAL, ALLOCATABLE :: CH4RVMR_tmp(:, :)
61  REAL CH4RVMR_glo(klon_glo,klev)
62!
63  INTEGER i, k, kk, j
64 
65! For NetCDF:
66  INTEGER ncid_in  ! IDs for input files
67  INTEGER varid, ncerr
68   
69  INTEGER, PARAMETER :: lev_input=17
70!--pressure at interfaces of input data (in Pa)
71  REAL, DIMENSION(lev_input+1), PARAMETER ::          &
72                    paprs_input=(/                    &
73  1.00000002e+05,   6.06530673e+04,   3.67879449e+04, &
74  2.23130165e+04,   1.35335286e+04,   8.20850004e+03, &
75  4.97870695e+03,   3.01973841e+03,   1.83156393e+03, &
76  1.11089968e+03,   6.73794715e+02,   4.08677153e+02, &
77  2.47875223e+02,   1.50343923e+02,   9.11881985e+01, &
78  5.53084382e+01,   3.35462635e+01,   0.0           /)
79!
80 
81  IF (debutphy .OR. mth_cur .NE. mth_pre) THEN
82     
83!--preparation of global fields
84     CALL gather(paprs, paprs_glo)
85     
86     IF (is_mpi_root.AND.is_omp_root) THEN
87       
88        print *,'In stratH2O_methox read file: ch4r_annual_lmdz.nc'
89       
90        !--reading strat. H2O source (CH4 oxidation) files
91        CALL nf95_open("ch4r_annual_lmdz.nc", nf90_nowrite, ncid_in)
92       
93        CALL  nf95_inq_varid(ncid_in, "LEV", varid)
94        CALL nf95_gw_var(ncid_in, varid, lev)
95        n_lev = size(lev)
96       
97        CALL nf95_inq_varid(ncid_in, "LAT", varid)
98        CALL nf95_gw_var(ncid_in, varid, latitude)
99        n_lat = size(latitude)
100       
101        CALL nf95_inq_varid(ncid_in, "LON", varid)
102        CALL nf95_gw_var(ncid_in, varid, longitude)
103        n_lon = size(longitude)
104       
105        CALL nf95_inq_varid(ncid_in, "TIME", varid)
106        CALL nf95_gw_var(ncid_in, varid, time)
107        n_mth = size(time)
108       
109        IF (.NOT.ALLOCATED(CH4RVMR_in)) ALLOCATE(CH4RVMR_in(n_lon, n_lat, n_lev, n_mth))
110       
111        CALL nf95_inq_varid(ncid_in, "CH4RVMR", varid)
112        ncerr = nf90_get_var(ncid_in, varid, CH4RVMR_in)
113        print *,'code erreur CH4RVMR=', ncerr, varid
114       
115        CALL nf95_close(ncid_in)
116       
117        IF (.NOT.ALLOCATED(CH4RVMR_mth)) ALLOCATE(CH4RVMR_mth(n_lon, n_lat, n_lev))
118        IF (.NOT.ALLOCATED(CH4RVMR_tmp)) ALLOCATE(CH4RVMR_tmp(klon_glo, n_lev))
119       
120        !---select the correct month,
121        !---correct latitudinal order,convert input from volume mixing ratio to mass mixing ratio
122        DO j=1,n_lat
123           ! convert VMR/s in MMR/s (mass mixing ratio/s or kg H2O/kg air/s)
124           ! x2 because CH4->2*H2O
125           CH4RVMR_mth(:,j,:) = 2*CH4RVMR_in(:,n_lat+1-j,:,mth_cur)*mH2Omol/mAIRmol
126        ENDDO
127       
128        !---reduce to a klon_glo grid but keep the levels
129        CALL grid2dTo1d_glo(CH4RVMR_mth,CH4RVMR_tmp)
130       
131        !---regrid weighted climatologies
132        DO i=1, klon_glo
133           DO k=1, klev
134             
135              CH4RVMR_glo(i,k)=0.0
136             
137              DO kk=1, n_lev
138                 !
139                 CH4RVMR_glo(i,k)=CH4RVMR_glo(i,k)+ &
140                      MAX(0.0,MIN(paprs_glo(i,k),paprs_input(kk)) &
141                      -MAX(paprs_glo(i,k+1),paprs_input(kk+1))) &
142                      *CH4RVMR_tmp(i,kk)/(paprs_glo(i,k)-paprs_glo(i,k+1))
143              ENDDO ! kk loop
144           ENDDO ! k loop
145        ENDDO ! i loop
146     ENDIF !--is_mpi_root and is_omp_root
147     
148     !--keep memory of previous month
149     mth_pre=mth_cur
150     
151     !--scatter global fields around
152     CALL scatter(CH4RVMR_glo, dq_ch4mmr)
153     
154     IF (is_mpi_root.AND.is_omp_root) THEN
155        DEALLOCATE(CH4RVMR_in,CH4RVMR_mth,CH4RVMR_tmp)
156     ENDIF !--is_mpi_root and is_omp_root
157  ENDIF ! debutphy.OR.new month
158 
159  RETURN
160 
161END SUBROUTINE stratH2O_methox
Note: See TracBrowser for help on using the repository browser.