source: LMDZ6/branches/LMDZ-QUEST/libf/phylmd/rrtm/readaerosolstrato1_rrtm.F90 @ 5429

Last change on this file since 5429 was 2744, checked in by oboucher, 8 years ago

Adding threadprivate for mth_pre variable
This should fix some reproducibility effects for runs longer than 1 months....

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