source: dynamico_lmdz/aquaplanet/LMDZ5/libf/phylmd/readaerosolstrato.F90 @ 3983

Last change on this file since 3983 was 3819, checked in by ymipsl, 10 years ago

Removed all iim et jjm depedency. Replaced by nbp_lon and nbp_lat.
Supress gr_fi_ecrit, replaced by grid1dTo2d_glo

YM

File size: 6.5 KB
Line 
1subroutine readaerosolstrato(debut)
2
3    use netcdf95, only: nf95_close, nf95_gw_var, nf95_inq_dimid, &
4                        nf95_inq_varid, nf95_open
5    use netcdf, only: nf90_get_var, nf90_noerr, nf90_nowrite
6
7    USE phys_cal_mod, ONLY : mth_cur
8    USE mod_grid_phy_lmdz, ONLY : nbp_lon, nbp_lat, klon_glo, grid2dTo1d_glo
9    USE mod_phys_lmdz_mpi_data, ONLY :  is_mpi_root
10    USE mod_phys_lmdz_para
11    USE phys_state_var_mod
12    USE phys_local_var_mod
13    USE aero_mod
14    USE dimphy
15
16    implicit none
17
18    include "YOMCST.h"
19
20! Variable input
21    logical debut
22
23! Variables locales
24    integer n_lat   ! number of latitudes in the input data
25    integer n_lon   ! number of longitudes in the input data
26    integer n_lev   ! number of levels in the input data
27    integer n_month ! number of months in the input data
28    real, pointer:: latitude(:)
29    real, pointer:: longitude(:)
30    real, pointer:: time(:)
31    real, pointer:: lev(:)
32    integer i, k, band, wave
33    integer, save :: mth_pre
34!$OMP THREADPRIVATE(mth_pre)
35
36    real, allocatable, dimension(:,:), save :: tau_aer_strat
37!$OMP THREADPRIVATE(tau_aer_strat)
38
39! Champs reconstitues
40    real, allocatable:: tauaerstrat(:, :, :, :)
41    real, allocatable:: tauaerstrat_mois(:, :, :)
42    real, allocatable:: tauaerstrat_mois_glo(:, :)
43
44! For NetCDF:
45    integer ncid_in  ! IDs for input files
46    integer varid, ncerr
47
48! Stratospheric aerosols optical properties
49! alpha_strat over the 2 bands is normalised by the 550 nm extinction coefficient
50! alpha_strat_wave is *not* normalised by the 550 nm extinction coefficient
51    real, dimension(nbands) :: alpha_strat, piz_strat, cg_strat
52    data alpha_strat/0.9922547, 0.7114912 /
53    data piz_strat  /0.9999998, 0.99762493/
54    data cg_strat   /0.73107845,0.73229635/
55    real, dimension(nwave) :: alpha_strat_wave
56    data alpha_strat_wave/3.36780953,3.34667683,3.20444202,3.0293026,2.82108808/
57
58!--------------------------------------------------------
59
60    IF (.not.ALLOCATED(tau_aer_strat)) ALLOCATE(tau_aer_strat(klon,klev))
61
62!--only read file if beginning of run or start of new month
63    IF (debut.OR.mth_cur.NE.mth_pre) THEN
64
65    IF (is_mpi_root) THEN
66
67    IF (nbands.NE.2) THEN
68        print *,'nbands doit etre egal a 2 dans readaerosolstrat'
69        STOP
70    ENDIF
71
72    CALL nf95_open("taustrat.nc", nf90_nowrite, ncid_in)
73
74    CALL nf95_inq_varid(ncid_in, "LEV", varid)
75    CALL nf95_gw_var(ncid_in, varid, lev)
76    n_lev = size(lev)
77    IF (n_lev.NE.klev) THEN
78       print *,'Le nombre de niveaux n est pas egal a klev'
79       STOP
80    ENDIF
81
82    CALL nf95_inq_varid(ncid_in, "LAT", varid)
83    CALL nf95_gw_var(ncid_in, varid, latitude)
84    n_lat = size(latitude)
85    print *, 'LAT aerosol strato=', n_lat, latitude
86    IF (n_lat.NE.nbp_lat) THEN
87       print *,'Le nombre de lat n est pas egal a nbp_lat'
88       STOP
89    ENDIF
90
91    CALL nf95_inq_varid(ncid_in, "LON", varid)
92    CALL nf95_gw_var(ncid_in, varid, longitude)
93    n_lon = size(longitude)
94    print *, 'LON aerosol strato=', n_lon, longitude
95    IF (n_lon.NE.nbp_lon) THEN
96       print *,'Le nombre de lon n est pas egal a nbp_lon'
97       STOP
98    ENDIF
99
100    CALL nf95_inq_varid(ncid_in, "TIME", varid)
101    CALL nf95_gw_var(ncid_in, varid, time)
102    n_month = size(time)
103    print *, 'TIME aerosol strato=', n_month, time
104    IF (n_month.NE.12) THEN
105       print *,'Le nombre de month n est pas egal a 12'
106       STOP
107    ENDIF
108
109    IF (.not.ALLOCATED(tauaerstrat))          ALLOCATE(tauaerstrat(n_lon, n_lat, n_lev, n_month))
110    IF (.not.ALLOCATED(tauaerstrat_mois))     ALLOCATE(tauaerstrat_mois(n_lon, n_lat, n_lev))
111    IF (.not.ALLOCATED(tauaerstrat_mois_glo)) ALLOCATE(tauaerstrat_mois_glo(klon_glo, n_lev))
112
113!--reading stratospheric AOD at 550 nm
114    CALL nf95_inq_varid(ncid_in, "TAUSTRAT", varid)
115    ncerr = nf90_get_var(ncid_in, varid, tauaerstrat)
116    print *,'code erreur readaerosolstrato=', ncerr, varid
117
118    CALL nf95_close(ncid_in)
119
120!---select the correct month
121    IF (mth_cur.LT.1.OR.mth_cur.GT.12) THEN
122      print *,'probleme avec le mois dans readaerosolstrat =', mth_cur
123    ENDIF
124    tauaerstrat_mois(:,:,:) = tauaerstrat(:,:,:,mth_cur)
125
126!---reduce to a klon_glo grid
127    CALL grid2dTo1d_glo(tauaerstrat_mois,tauaerstrat_mois_glo)
128
129    ENDIF !--is_mpi_root
130
131!--scatter on all proc
132    CALL scatter(tauaerstrat_mois_glo,tau_aer_strat)
133
134!--keep memory of previous month
135    mth_pre=mth_cur
136!
137    IF (is_mpi_root) THEN
138!
139    DEALLOCATE(tauaerstrat)
140    DEALLOCATE(tauaerstrat_mois)
141    DEALLOCATE(tauaerstrat_mois_glo)
142!
143    ENDIF !-is_mpi_root
144
145    ENDIF !--debut ou nouveau mois
146
147!--total vertical aod at the 6 wavelengths
148    DO wave=1, nwave
149    DO k=1, klev
150    tausum_aero(:,wave,id_STRAT_phy)=tausum_aero(:,wave,id_STRAT_phy)+tau_aer_strat(:,k)*alpha_strat_wave(wave)/alpha_strat_wave(2)
151    ENDDO
152    ENDDO
153
154!--weighted average for cg, piz and tau, adding strat aerosols on top of tropospheric ones
155    DO band=1, nbands
156!--anthropogenic aerosols bands 1 and 2
157    cg_aero(:,:,3,band)  = ( cg_aero(:,:,3,band)*piz_aero(:,:,3,band)*tau_aero(:,:,3,band) +          &
158                             cg_strat(band)*piz_strat(band)*alpha_strat(band)*tau_aer_strat(:,:) ) /  &
159                             MAX( piz_aero(:,:,3,band)*tau_aero(:,:,3,band) +                         &
160                                  piz_strat(band)*alpha_strat(band)*tau_aer_strat(:,:), 1.e-15 )
161    piz_aero(:,:,3,band)  = ( piz_aero(:,:,3,band)*tau_aero(:,:,3,band) +                             &
162                              piz_strat(band)*alpha_strat(band)*tau_aer_strat(:,:) ) /                &
163                              MAX( tau_aero(:,:,3,band) + alpha_strat(band)*tau_aer_strat(:,:), 1.e-15 )
164    tau_aero(:,:,3,band)  = tau_aero(:,:,3,band) + alpha_strat(band)*tau_aer_strat(:,:)
165!--natural aerosols bands 1 and 2
166    cg_aero(:,:,2,band)  = ( cg_aero(:,:,2,band)*piz_aero(:,:,2,band)*tau_aero(:,:,2,band) +          &
167                             cg_strat(band)*piz_strat(band)*alpha_strat(band)*tau_aer_strat(:,:) ) /  &
168                             MAX( piz_aero(:,:,2,band)*tau_aero(:,:,2,band) +                         &
169                                  piz_strat(band)*alpha_strat(band)*tau_aer_strat(:,:), 1.e-15 )
170    piz_aero(:,:,2,band)  = ( piz_aero(:,:,2,band)*tau_aero(:,:,2,band) +                             &
171                              piz_strat(band)*alpha_strat(band)*tau_aer_strat(:,:) ) /                &
172                              MAX( tau_aero(:,:,2,band) + alpha_strat(band)*tau_aer_strat(:,:),1.e-15 )
173    tau_aero(:,:,2,band)  = tau_aero(:,:,2,band) + alpha_strat(band)*tau_aer_strat(:,:)
174    ENDDO
175
176end subroutine readaerosolstrato
Note: See TracBrowser for help on using the repository browser.