source: lmdz_wrf/WRFV3/lmdz/readaerosolstrato.F90 @ 1

Last change on this file since 1 was 1, checked in by lfita, 10 years ago
  • -- --- Opening of the WRF+LMDZ coupling repository --- -- -

WRF: version v3.3
LMDZ: version v1818

More details in:

File size: 6.3 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
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    real, allocatable:: tauaerstrat_mois_glo_bands(:,:,:)
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    IF (is_mpi_root) THEN
64
65    IF (debut.OR.mth_cur.NE.mth_pre) 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.jjm+1) THEN
87       print *,'Le nombre de lat n est pas egal a jjm+1'
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.iim) THEN
96       print *,'Le nombre de lon n est pas egal a iim'
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    ALLOCATE(tauaerstrat(n_lon, n_lat, n_lev, n_month))
110    ALLOCATE(tauaerstrat_mois(n_lon, n_lat, n_lev))
111    ALLOCATE(tauaerstrat_mois_glo(klon_glo, n_lev))
112    ALLOCATE(tauaerstrat_mois_glo_bands(klon_glo, n_lev,nbands))
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!--scatter on all proc
131    CALL scatter(tauaerstrat_mois_glo,tau_aer_strat)
132
133    DEALLOCATE(tauaerstrat)
134    DEALLOCATE(tauaerstrat_mois)
135    DEALLOCATE(tauaerstrat_mois_glo)
136 
137    mth_pre=mth_cur
138
139    ENDIF !--debut ou nouveau mois
140
141    ENDIF !--is_mpi_root
142
143!--total vertical aod at the 5 wavelengths
144    DO wave=1, nwave
145    DO k=1, klev
146    tausum_aero(:,wave,id_strat)=tausum_aero(:,wave,id_strat)+tau_aer_strat(:,k)*alpha_strat_wave(wave)/alpha_strat_wave(2)
147    ENDDO
148    ENDDO
149
150!--weighted average for cg, piz and tau, adding strat aerosols on top of tropospheric ones
151    DO band=1, nbands
152!--anthropogenic aerosols bands 1 and 2
153    cg_aero(:,:,3,band)  = ( cg_aero(:,:,3,band)*piz_aero(:,:,3,band)*tau_aero(:,:,3,band) +          &
154                             cg_strat(band)*piz_strat(band)*alpha_strat(band)*tau_aer_strat(:,:) ) /  &
155                             MAX( piz_aero(:,:,3,band)*tau_aero(:,:,3,band) +                         &
156                                  piz_strat(band)*alpha_strat(band)*tau_aer_strat(:,:), 1.e-15 )
157    piz_aero(:,:,3,band)  = ( piz_aero(:,:,3,band)*tau_aero(:,:,3,band) +                             &
158                              piz_strat(band)*alpha_strat(band)*tau_aer_strat(:,:) ) /                &
159                              MAX( tau_aero(:,:,3,band) + alpha_strat(band)*tau_aer_strat(:,:), 1.e-15 )
160    tau_aero(:,:,3,band)  = tau_aero(:,:,3,band) + alpha_strat(band)*tau_aer_strat(:,:)
161!--natural aerosols bands 1 and 2
162    cg_aero(:,:,2,band)  = ( cg_aero(:,:,2,band)*piz_aero(:,:,2,band)*tau_aero(:,:,2,band) +          &
163                             cg_strat(band)*piz_strat(band)*alpha_strat(band)*tau_aer_strat(:,:) ) /  &
164                             MAX( piz_aero(:,:,2,band)*tau_aero(:,:,2,band) +                         &
165                                  piz_strat(band)*alpha_strat(band)*tau_aer_strat(:,:), 1.e-15 )
166    piz_aero(:,:,2,band)  = ( piz_aero(:,:,2,band)*tau_aero(:,:,2,band) +                             &
167                              piz_strat(band)*alpha_strat(band)*tau_aer_strat(:,:) ) /                &
168                              MAX( tau_aero(:,:,2,band) + alpha_strat(band)*tau_aer_strat(:,:),1.e-15 )
169    tau_aero(:,:,2,band)  = tau_aero(:,:,2,band) + alpha_strat(band)*tau_aer_strat(:,:)
170    ENDDO
171
172end subroutine readaerosolstrato
Note: See TracBrowser for help on using the repository browser.