source: LMDZ5/trunk/libf/phylmd/readaerosolstrato.F90 @ 2306

Last change on this file since 2306 was 2152, checked in by fhourdin, 10 years ago

Corrections cosmétiques pour RRTM
Bug fixing for RRTM

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