source: LMDZ5/branches/testing/libf/phylmd/rrtm/readaerosolstrato_rrtm.F90 @ 2408

Last change on this file since 2408 was 2408, checked in by Laurent Fairhead, 9 years ago

Merged trunk changes r2298:2396 into testing branch

  • Property copyright set to
    Name of program: LMDZ
    Creation date: 1984
    Version: LMDZ5
    License: CeCILL version 2
    Holder: Laboratoire de m\'et\'eorologie dynamique, CNRS, UMR 8539
    See the license file in the root directory
  • Property svn:keywords set to Author Date Id Revi
File size: 9.1 KB
Line 
1!
2! $Id: readaerosolstrato_rrtm.F90 2408 2015-12-14 10:43:09Z fairhead $
3!
4subroutine readaerosolstrato_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_para
14    USE phys_state_var_mod
15    USE phys_local_var_mod
16    USE aero_mod
17    USE dimphy
18    USE YOERAD   , ONLY : NLW
19
20    implicit none
21
22    include "YOMCST.h"
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) :: 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
66    real :: alpha_lw_strat_wave
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    IF (debut.OR.mth_cur.NE.mth_pre) THEN
82
83    IF (is_mpi_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
148
149!--keep memory of previous month
150    mth_pre=mth_cur
151
152!--scatter on all proc
153    CALL scatter(tauaerstrat_mois_glo,tau_aer_strat)
154
155    IF (is_mpi_root) THEN
156!
157    DEALLOCATE(tauaerstrat)
158    DEALLOCATE(tauaerstrat_mois)
159    DEALLOCATE(tauaerstrat_mois_glo)
160!
161    ENDIF !--is_mpi_root
162
163    ENDIF !--debut ou nouveau mois
164
165!--total vertical aod at the 5 SW wavelengths
166    DO wave=1, nwave
167    DO k=1, klev
168    tausum_aero(:,wave,id_STRAT_phy)=tausum_aero(:,wave,id_STRAT_phy)+ &
169       tau_aer_strat(:,k)*alpha_sw_strat_wave(wave)/alpha_sw_strat_wave(2)
170    ENDDO
171    ENDDO
172
173!--weighted average for cg, piz and tau, adding strat aerosols on top of tropospheric ones
174    DO band=1, nbands_sw_rrtm
175!--anthropogenic aerosols bands 1 to nbands_sw_rrtm
176    cg_aero_sw_rrtm(:,:,2,band)  = ( cg_aero_sw_rrtm(:,:,2,band)*piz_aero_sw_rrtm(:,:,2,band)*tau_aero_sw_rrtm(:,:,2,band) + &
177                                  cg_sw_strat(band)*piz_sw_strat(band)*alpha_sw_strat(band)*tau_aer_strat(:,:) ) /           &
178                             MAX( piz_aero_sw_rrtm(:,:,2,band)*tau_aero_sw_rrtm(:,:,2,band) +                                &
179                                  piz_sw_strat(band)*alpha_sw_strat(band)*tau_aer_strat(:,:), 1.e-15 )
180    piz_aero_sw_rrtm(:,:,2,band)  = ( piz_aero_sw_rrtm(:,:,2,band)*tau_aero_sw_rrtm(:,:,2,band) +                            &
181                              piz_sw_strat(band)*alpha_sw_strat(band)*tau_aer_strat(:,:) ) /                                 &
182                              MAX( tau_aero_sw_rrtm(:,:,2,band) + alpha_sw_strat(band)*tau_aer_strat(:,:), 1.e-15 )
183    tau_aero_sw_rrtm(:,:,2,band)  = tau_aero_sw_rrtm(:,:,2,band) + alpha_sw_strat(band)*tau_aer_strat(:,:)
184!--natural aerosols bands 1 to nbands_sw_rrtm
185    cg_aero_sw_rrtm(:,:,1,band)  = ( cg_aero_sw_rrtm(:,:,1,band)*piz_aero_sw_rrtm(:,:,1,band)*tau_aero_sw_rrtm(:,:,1,band) + &
186                             cg_sw_strat(band)*piz_sw_strat(band)*alpha_sw_strat(band)*tau_aer_strat(:,:) ) /                &
187                             MAX( piz_aero_sw_rrtm(:,:,1,band)*tau_aero_sw_rrtm(:,:,1,band) +                                &
188                                  piz_sw_strat(band)*alpha_sw_strat(band)*tau_aer_strat(:,:), 1.e-15 )
189    piz_aero_sw_rrtm(:,:,1,band)  = ( piz_aero_sw_rrtm(:,:,1,band)*tau_aero_sw_rrtm(:,:,1,band) +                            &
190                              piz_sw_strat(band)*alpha_sw_strat(band)*tau_aer_strat(:,:) ) /                                 &
191                              MAX( tau_aero_sw_rrtm(:,:,1,band) + alpha_sw_strat(band)*tau_aer_strat(:,:),1.e-15 )
192    tau_aero_sw_rrtm(:,:,1,band)  = tau_aero_sw_rrtm(:,:,1,band) + alpha_sw_strat(band)*tau_aer_strat(:,:)
193!--no stratospheric aerosol in index 1 for these tests
194!    cg_aero_sw_rrtm(:,:,1,band)  =  cg_aero_sw_rrtm(:,:,1,band)
195!    piz_aero_sw_rrtm(:,:,1,band)  = piz_aero_sw_rrtm(:,:,1,band)
196!    tau_aero_sw_rrtm(:,:,1,band)  = tau_aero_sw_rrtm(:,:,1,band)
197    ENDDO
198
199!--stratospheric AOD in LW
200    IF (nbands_lw_rrtm .NE. NLW) then
201      print*, 'different values for NLW (=',NLW,') and nbands_lw_rrtm (=', nbands_lw_rrtm, ')'
202      STOP
203    ENDIF
204
205    DO band=1, nbands_lw_rrtm
206    tau_aero_lw_rrtm(:,:,2,band)  = tau_aero_lw_rrtm(:,:,2,band) + alpha_lw_abs_rrtm(band)*tau_aer_strat(:,:)
207    tau_aero_lw_rrtm(:,:,1,band)  = tau_aero_lw_rrtm(:,:,1,band) + alpha_lw_abs_rrtm(band)*tau_aer_strat(:,:)
208!--no stratospheric aerosols in index 1 for these tests
209!    tau_aero_lw_rrtm(:,:,1,band)  = tau_aero_lw_rrtm(:,:,1,band)
210    ENDDO
211
212!--default SSA value if there is no aerosol
213!--to avoid 0 values that seems to cause some problem to RRTM
214    WHERE (tau_aero_sw_rrtm.LT.1.e-14)
215      piz_aero_sw_rrtm = 1.0
216    ENDWHERE
217
218!--in principle this should not be necessary
219!--as these variables have min values already but just in case
220!--put 1e-15 min value to both SW and LW AOD
221    tau_aero_sw_rrtm = MAX(tau_aero_sw_rrtm,1.e-15)
222    tau_aero_lw_rrtm = MAX(tau_aero_lw_rrtm,1.e-15)
223
224end subroutine readaerosolstrato_rrtm
Note: See TracBrowser for help on using the repository browser.