source: LMDZ5/trunk/libf/phylmd/rrtm/readaerosolstrato1_rrtm.F90 @ 2724

Last change on this file since 2724 was 2694, checked in by oboucher, 8 years ago

Cosmetic changes

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