source: LMDZ5/branches/IPSLCM5A2.1_ISO/libf/phyiso/rrtm/readaerosolstrato1_rrtm.F90 @ 4249

Last change on this file since 4249 was 3331, checked in by acozic, 7 years ago

Add modification for isotopes

  • Property svn:executable set to *
File size: 9.7 KB
Line 
1!
2! $Id: readaerosolstrato1_rrtm.F90 2526 2016-05-26 22:13:40Z oboucher $
3!
4subroutine readaerosolstrato1_rrtm(debut)
5
6    use netcdf95, only: nf95_close, nf95_gw_var, nf95_inq_dimid, &
7                        nf95_inq_varid, nf95_open
8    use netcdf, only: nf90_get_var, nf90_noerr, nf90_nowrite
9
10    USE phys_cal_mod, ONLY : mth_cur
11    USE mod_grid_phy_lmdz, ONLY: nbp_lon, nbp_lat, klon_glo, grid2dTo1d_glo
12    USE mod_phys_lmdz_mpi_data, ONLY :  is_mpi_root
13    USE mod_phys_lmdz_omp_data, ONLY :  is_omp_root
14    USE mod_phys_lmdz_para
15    USE phys_state_var_mod
16    USE phys_local_var_mod
17    USE aero_mod
18    USE dimphy
19    USE YOERAD   , ONLY : NLW
20
21    implicit none
22
23    include "YOMCST.h"
24
25! Variable input
26    logical debut
27
28! Variables locales
29    integer n_lat   ! number of latitudes in the input data
30    integer n_lon   ! number of longitudes in the input data
31    integer n_lev   ! number of levels in the input data
32    integer n_month ! number of months in the input data
33    real, pointer:: latitude(:)
34    real, pointer:: longitude(:)
35    real, pointer:: time(:)
36    real, pointer:: lev(:)
37    integer k, band, wave, i
38    integer, save :: mth_pre
39
40    real, allocatable, dimension(:,:), save :: tau_aer_strat
41!$OMP THREADPRIVATE(tau_aer_strat)
42
43! Champs reconstitues
44    real, allocatable:: tauaerstrat(:, :, :, :)
45    real, allocatable:: tauaerstrat_mois(:, :, :)
46    real, allocatable:: tauaerstrat_mois_glo(:, :)
47
48    real, allocatable:: sum_tau_aer_strat(:)
49
50! For NetCDF:
51    integer ncid_in  ! IDs for input files
52    integer varid, ncerr
53
54! Stratospheric aerosols optical properties
55! alpha_sw_strat over the 6 bands is normalised by the 550 nm extinction coefficient
56    real, dimension(nbands_sw_rrtm) :: alpha_sw_strat, piz_sw_strat, cg_sw_strat
57    data alpha_sw_strat/0.8545564, 0.8451642, 0.9821724, 0.8145110, 0.3073565, 7.7966176E-02/
58    data cg_sw_strat   /0.6997170, 0.6810035, 0.7403592, 0.7562674, 0.6676504, 0.3478689/
59    data piz_sw_strat  /0.9999998, 0.9999998, 1.000000000, 0.9999958, 0.9977155, 0.4510679/
60!
61!--diagnostics AOD in the SW
62! alpha_sw_strat_wave is *not* normalised by the 550 nm extinction coefficient
63    real, dimension(nwave_sw) :: alpha_sw_strat_wave
64    data alpha_sw_strat_wave/3.708007,4.125824,4.136584,3.887478,3.507738/
65!
66!--diagnostics AOD in the LW at 10 um (not normalised by the 550 nm ext coefficient
67    real :: alpha_lw_strat_wave(nwave_lw)
68    data alpha_lw_strat_wave/0.2746812/
69!
70    real, dimension(nbands_lw_rrtm) :: alpha_lw_abs_rrtm
71    data alpha_lw_abs_rrtm/   8.8340312E-02, 6.9856711E-02, 6.2652975E-02, 5.7188231E-02, &
72                              6.3157059E-02, 5.5072524E-02, 5.0571125E-02, 0.1349073, &   
73                              0.1381676, 9.6506312E-02, 5.1312990E-02, 2.4256418E-02, &
74                              2.7191756E-02, 3.3862915E-02, 1.6132960E-02, 1.4275438E-02/ ! calculated with Mie_SW_LW_RRTM_V2.4 (bimodal, corrected)
75                                                                                          ! for r_0=/0.13E-6, 0.41E-6/ m, sigma_g=/1.26, 1.30/
76                                                                                          ! order: increasing wavelength!
77!--------------------------------------------------------
78
79    IF (.not.ALLOCATED(tau_aer_strat)) ALLOCATE(tau_aer_strat(klon,klev))
80    IF (.not.ALLOCATED(sum_tau_aer_strat)) ALLOCATE(sum_tau_aer_strat(klon))
81
82!--we only read monthly strat aerosol data
83    IF (debut.OR.mth_cur.NE.mth_pre) THEN
84
85!--only root reads the data
86    IF (is_mpi_root.AND.is_omp_root) THEN
87
88    IF (nbands_sw_rrtm.NE.6) THEN
89        print *,'nbands_sw_rrtm doit etre egal a 6 dans readaerosolstrat_rrtm'
90        STOP
91    ENDIF
92
93    CALL nf95_open("taustrat.nc", nf90_nowrite, ncid_in)
94
95    CALL nf95_inq_varid(ncid_in, "LEV", varid)
96    CALL nf95_gw_var(ncid_in, varid, lev)
97    n_lev = size(lev)
98    IF (n_lev.NE.klev) THEN
99       print *,'Le nombre de niveaux n est pas egal a klev'
100       STOP
101    ENDIF
102
103    CALL nf95_inq_varid(ncid_in, "LAT", varid)
104    CALL nf95_gw_var(ncid_in, varid, latitude)
105    n_lat = size(latitude)
106    print *, 'LAT aerosol strato=', n_lat, latitude
107    IF (n_lat.NE.nbp_lat) THEN
108       print *,'Le nombre de lat n est pas egal a nbp_lat'
109       STOP
110    ENDIF
111
112    CALL nf95_inq_varid(ncid_in, "LON", varid)
113    CALL nf95_gw_var(ncid_in, varid, longitude)
114    n_lon = size(longitude)
115    print *, 'LON aerosol strato=', n_lon, longitude
116    IF (n_lon.NE.nbp_lon) THEN
117       print *,'Le nombre de lon n est pas egal a nbp_lon'
118       STOP
119    ENDIF
120
121    CALL nf95_inq_varid(ncid_in, "TIME", varid)
122    CALL nf95_gw_var(ncid_in, varid, time)
123    n_month = size(time)
124    print *, 'TIME aerosol strato=', n_month, time
125    IF (n_month.NE.12) THEN
126       print *,'Le nombre de month n est pas egal a 12'
127       STOP
128    ENDIF
129
130    ALLOCATE(tauaerstrat(n_lon, n_lat, n_lev, n_month))
131    ALLOCATE(tauaerstrat_mois(n_lon, n_lat, n_lev))
132    ALLOCATE(tauaerstrat_mois_glo(klon_glo, n_lev))
133
134!--reading stratospheric AOD at 550 nm
135    CALL nf95_inq_varid(ncid_in, "TAUSTRAT", varid)
136    ncerr = nf90_get_var(ncid_in, varid, tauaerstrat)
137    print *,'code erreur readaerosolstrato=', ncerr, varid
138
139    CALL nf95_close(ncid_in)
140
141!---select the correct month
142    IF (mth_cur.LT.1.OR.mth_cur.GT.12) THEN
143      print *,'probleme avec le mois dans readaerosolstrat =', mth_cur
144    ENDIF
145    tauaerstrat_mois(:,:,:) = tauaerstrat(:,:,:,mth_cur)
146
147!---reduce to a klon_glo grid
148    CALL grid2dTo1d_glo(tauaerstrat_mois,tauaerstrat_mois_glo)
149
150    ENDIF !--is_mpi_root and is_omp_root
151
152!$OMP BARRIER
153
154!--keep memory of previous month
155    mth_pre=mth_cur
156
157!--scatter on all proc
158    CALL scatter(tauaerstrat_mois_glo,tau_aer_strat)
159
160    IF (is_mpi_root.AND.is_omp_root) THEN
161!
162    DEALLOCATE(tauaerstrat)
163    DEALLOCATE(tauaerstrat_mois)
164    DEALLOCATE(tauaerstrat_mois_glo)
165!
166    ENDIF !--is_mpi_root and is_omp_root
167
168!$OMP BARRIER
169
170    ENDIF !--debut ou nouveau mois
171
172!--total vertical aod at the 5 SW wavelengths
173    DO wave=1, nwave_sw
174    DO k=1, klev
175      tausum_aero(:,wave,id_STRAT_phy)=tausum_aero(:,wave,id_STRAT_phy)+ &
176          tau_aer_strat(:,k)*alpha_sw_strat_wave(wave)/alpha_sw_strat_wave(2)
177    ENDDO
178    ENDDO
179
180!--weighted average for cg, piz and tau, adding strat aerosols on top of tropospheric ones
181    DO band=1, nbands_sw_rrtm
182!--anthropogenic aerosols bands 1 to nbands_sw_rrtm
183    cg_aero_sw_rrtm(:,:,2,band)  = ( cg_aero_sw_rrtm(:,:,2,band)*piz_aero_sw_rrtm(:,:,2,band)*tau_aero_sw_rrtm(:,:,2,band) + &
184                                  cg_sw_strat(band)*piz_sw_strat(band)*alpha_sw_strat(band)*tau_aer_strat(:,:) ) /           &
185                             MAX( piz_aero_sw_rrtm(:,:,2,band)*tau_aero_sw_rrtm(:,:,2,band) +                                &
186                                  piz_sw_strat(band)*alpha_sw_strat(band)*tau_aer_strat(:,:), 1.e-15 )
187    piz_aero_sw_rrtm(:,:,2,band)  = ( piz_aero_sw_rrtm(:,:,2,band)*tau_aero_sw_rrtm(:,:,2,band) +                            &
188                              piz_sw_strat(band)*alpha_sw_strat(band)*tau_aer_strat(:,:) ) /                                 &
189                              MAX( tau_aero_sw_rrtm(:,:,2,band) + alpha_sw_strat(band)*tau_aer_strat(:,:), 1.e-15 )
190    tau_aero_sw_rrtm(:,:,2,band)  = tau_aero_sw_rrtm(:,:,2,band) + alpha_sw_strat(band)*tau_aer_strat(:,:)
191!--natural aerosols bands 1 to nbands_sw_rrtm
192    cg_aero_sw_rrtm(:,:,1,band)  = ( cg_aero_sw_rrtm(:,:,1,band)*piz_aero_sw_rrtm(:,:,1,band)*tau_aero_sw_rrtm(:,:,1,band) + &
193                             cg_sw_strat(band)*piz_sw_strat(band)*alpha_sw_strat(band)*tau_aer_strat(:,:) ) /                &
194                             MAX( piz_aero_sw_rrtm(:,:,1,band)*tau_aero_sw_rrtm(:,:,1,band) +                                &
195                                  piz_sw_strat(band)*alpha_sw_strat(band)*tau_aer_strat(:,:), 1.e-15 )
196    piz_aero_sw_rrtm(:,:,1,band)  = ( piz_aero_sw_rrtm(:,:,1,band)*tau_aero_sw_rrtm(:,:,1,band) +                            &
197                              piz_sw_strat(band)*alpha_sw_strat(band)*tau_aer_strat(:,:) ) /                                 &
198                              MAX( tau_aero_sw_rrtm(:,:,1,band) + alpha_sw_strat(band)*tau_aer_strat(:,:),1.e-15 )
199    tau_aero_sw_rrtm(:,:,1,band)  = tau_aero_sw_rrtm(:,:,1,band) + alpha_sw_strat(band)*tau_aer_strat(:,:)
200!--no stratospheric aerosol in index 1 for these tests
201!    cg_aero_sw_rrtm(:,:,1,band)  =  cg_aero_sw_rrtm(:,:,1,band)
202!    piz_aero_sw_rrtm(:,:,1,band)  = piz_aero_sw_rrtm(:,:,1,band)
203!    tau_aero_sw_rrtm(:,:,1,band)  = tau_aero_sw_rrtm(:,:,1,band)
204    ENDDO
205
206!--stratospheric AOD in LW
207    IF (nbands_lw_rrtm .NE. NLW) then
208      print*, 'different values for NLW (=',NLW,') and nbands_lw_rrtm (=', nbands_lw_rrtm, ')'
209      STOP
210    ENDIF
211
212!--total vertical aod at the 1 LW wavelength
213    DO wave=1, nwave_lw
214    DO k=1, klev
215      tausum_aero(:,nwave_sw+wave,id_STRAT_phy)=tausum_aero(:,nwave_sw+wave,id_STRAT_phy)+ &
216         tau_aer_strat(:,k)*alpha_lw_strat_wave(wave)/alpha_sw_strat_wave(2)
217    ENDDO
218    ENDDO
219
220    DO band=1, nbands_lw_rrtm
221    tau_aero_lw_rrtm(:,:,2,band)  = tau_aero_lw_rrtm(:,:,2,band) + alpha_lw_abs_rrtm(band)*tau_aer_strat(:,:)
222    tau_aero_lw_rrtm(:,:,1,band)  = tau_aero_lw_rrtm(:,:,1,band) + alpha_lw_abs_rrtm(band)*tau_aer_strat(:,:)
223!--no stratospheric aerosols in index 1 for these tests
224!    tau_aero_lw_rrtm(:,:,1,band)  = tau_aero_lw_rrtm(:,:,1,band)
225    ENDDO
226
227!--default SSA value if there is no aerosol
228!--to avoid 0 values that seems to cause some problem to RRTM
229    WHERE (tau_aero_sw_rrtm.LT.1.e-14)
230      piz_aero_sw_rrtm = 1.0
231    ENDWHERE
232
233!--in principle this should not be necessary
234!--as these variables have min values already but just in case
235!--put 1e-15 min value to both SW and LW AOD
236    tau_aero_sw_rrtm = MAX(tau_aero_sw_rrtm,1.e-15)
237    tau_aero_lw_rrtm = MAX(tau_aero_lw_rrtm,1.e-15)
238
239end subroutine readaerosolstrato1_rrtm
Note: See TracBrowser for help on using the repository browser.