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

Last change on this file since 2056 was 2056, checked in by Laurent Fairhead, 10 years ago

Merged trunk changes r1997:2055 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: 6.7 KB
Line 
1!
2! $Id: readaerosolstrato_rrtm.F90 2056 2014-06-11 13:46:46Z 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
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
19    implicit none
20
21    include "YOMCST.h"
22    include "dimensions.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
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    real, allocatable:: tauaerstrat_mois_glo_bands(:,:,:)
47
48! For NetCDF:
49    integer ncid_in  ! IDs for input files
50    integer varid, ncerr
51
52! Stratospheric aerosols optical properties
53! alpha_strat over the 2 bands is normalised by the 550 nm extinction coefficient
54! alpha_strat_wave is *not* normalised by the 550 nm extinction coefficient
55    real, dimension(nbands_rrtm) :: alpha_strat, piz_strat, cg_strat
56    data alpha_strat/0.938538969, 0.990073204, 0.992904723, 0.829215884, 0.439313501, 0.156857833/
57    data cg_strat   /0.699142992, 0.716326416, 0.735462785, 0.736726701, 0.712068975, 0.575097859/
58    data piz_strat  /1.000000000, 1.000000000, 1.000000000, 1.000000000, 0.997781098, 0.452584684/
59    real, dimension(nwave) :: alpha_strat_wave
60    data alpha_strat_wave/3.36780953,3.34667683,3.20444202,3.0293026,2.82108808/
61
62!--------------------------------------------------------
63
64    IF (.not.ALLOCATED(tau_aer_strat)) ALLOCATE(tau_aer_strat(klon,klev))
65
66    IF (is_mpi_root) THEN
67
68    IF (debut.OR.mth_cur.NE.mth_pre) THEN
69
70    IF (nbands_rrtm.NE.6) THEN
71        print *,'nbands_rrtm doit etre egal a 6 dans readaerosolstrat_rrtm'
72        STOP
73    ENDIF
74
75    CALL nf95_open("taustrat.nc", nf90_nowrite, ncid_in)
76
77    CALL nf95_inq_varid(ncid_in, "LEV", varid)
78    CALL nf95_gw_var(ncid_in, varid, lev)
79    n_lev = size(lev)
80    IF (n_lev.NE.klev) THEN
81       print *,'Le nombre de niveaux n est pas egal a klev'
82       STOP
83    ENDIF
84
85    CALL nf95_inq_varid(ncid_in, "LAT", varid)
86    CALL nf95_gw_var(ncid_in, varid, latitude)
87    n_lat = size(latitude)
88    print *, 'LAT aerosol strato=', n_lat, latitude
89    IF (n_lat.NE.jjm+1) THEN
90       print *,'Le nombre de lat n est pas egal a jjm+1'
91       STOP
92    ENDIF
93
94    CALL nf95_inq_varid(ncid_in, "LON", varid)
95    CALL nf95_gw_var(ncid_in, varid, longitude)
96    n_lon = size(longitude)
97    print *, 'LON aerosol strato=', n_lon, longitude
98    IF (n_lon.NE.iim) THEN
99       print *,'Le nombre de lon n est pas egal a iim'
100       STOP
101    ENDIF
102
103    CALL nf95_inq_varid(ncid_in, "TIME", varid)
104    CALL nf95_gw_var(ncid_in, varid, time)
105    n_month = size(time)
106    print *, 'TIME aerosol strato=', n_month, time
107    IF (n_month.NE.12) THEN
108       print *,'Le nombre de month n est pas egal a 12'
109       STOP
110    ENDIF
111
112    ALLOCATE(tauaerstrat(n_lon, n_lat, n_lev, n_month))
113    ALLOCATE(tauaerstrat_mois(n_lon, n_lat, n_lev))
114    ALLOCATE(tauaerstrat_mois_glo(klon_glo, n_lev))
115    ALLOCATE(tauaerstrat_mois_glo_bands(klon_glo, n_lev,nbands_rrtm))
116
117!--reading stratospheric AOD at 550 nm
118    CALL nf95_inq_varid(ncid_in, "TAUSTRAT", varid)
119    ncerr = nf90_get_var(ncid_in, varid, tauaerstrat)
120    print *,'code erreur readaerosolstrato=', ncerr, varid
121
122    CALL nf95_close(ncid_in)
123
124!---select the correct month
125    IF (mth_cur.LT.1.OR.mth_cur.GT.12) THEN
126      print *,'probleme avec le mois dans readaerosolstrat =', mth_cur
127    ENDIF
128    tauaerstrat_mois(:,:,:) = tauaerstrat(:,:,:,mth_cur)
129
130!---reduce to a klon_glo grid
131    CALL grid2dTo1d_glo(tauaerstrat_mois,tauaerstrat_mois_glo)
132
133!--scatter on all proc
134    CALL scatter(tauaerstrat_mois_glo,tau_aer_strat)
135
136    DEALLOCATE(tauaerstrat)
137    DEALLOCATE(tauaerstrat_mois)
138    DEALLOCATE(tauaerstrat_mois_glo)
139 
140    mth_pre=mth_cur
141
142    ENDIF !--debut ou nouveau mois
143
144    ENDIF !--is_mpi_root
145
146!--total vertical aod at the 5 wavelengths
147    DO wave=1, nwave
148    DO k=1, klev
149    tausum_aero(:,wave,id_strat)=tausum_aero(:,wave,id_strat)+tau_aer_strat(:,k)*alpha_strat_wave(wave)/alpha_strat_wave(2)
150    ENDDO
151    ENDDO
152
153!--weighted average for cg, piz and tau, adding strat aerosols on top of tropospheric ones
154    DO band=1, nbands_rrtm
155!--anthropogenic aerosols bands 1 to nbands_rrtm
156    cg_aero_rrtm(:,:,2,band)  = ( cg_aero_rrtm(:,:,2,band)*piz_aero_rrtm(:,:,2,band)*tau_aero_rrtm(:,:,2,band) + &
157                                  cg_strat(band)*piz_strat(band)*alpha_strat(band)*tau_aer_strat(:,:) ) /        &
158                             MAX( piz_aero_rrtm(:,:,2,band)*tau_aero_rrtm(:,:,2,band) +                          &
159                                  piz_strat(band)*alpha_strat(band)*tau_aer_strat(:,:), 1.e-15 )
160    piz_aero_rrtm(:,:,2,band)  = ( piz_aero_rrtm(:,:,2,band)*tau_aero_rrtm(:,:,2,band) +                         &
161                              piz_strat(band)*alpha_strat(band)*tau_aer_strat(:,:) ) /                           &
162                              MAX( tau_aero_rrtm(:,:,2,band) + alpha_strat(band)*tau_aer_strat(:,:), 1.e-15 )
163    tau_aero_rrtm(:,:,2,band)  = tau_aero_rrtm(:,:,2,band) + alpha_strat(band)*tau_aer_strat(:,:)
164!--natural aerosols bands 1 to nbands_rrtm
165    cg_aero_rrtm(:,:,1,band)  = ( cg_aero_rrtm(:,:,1,band)*piz_aero_rrtm(:,:,1,band)*tau_aero_rrtm(:,:,1,band) + &
166                             cg_strat(band)*piz_strat(band)*alpha_strat(band)*tau_aer_strat(:,:) ) /             &
167                             MAX( piz_aero_rrtm(:,:,1,band)*tau_aero_rrtm(:,:,1,band) +                          &
168                                  piz_strat(band)*alpha_strat(band)*tau_aer_strat(:,:), 1.e-15 )
169    piz_aero_rrtm(:,:,1,band)  = ( piz_aero_rrtm(:,:,1,band)*tau_aero_rrtm(:,:,1,band) +                         &
170                              piz_strat(band)*alpha_strat(band)*tau_aer_strat(:,:) ) /                           &
171                              MAX( tau_aero_rrtm(:,:,1,band) + alpha_strat(band)*tau_aer_strat(:,:),1.e-15 )
172    tau_aero_rrtm(:,:,1,band)  = tau_aero_rrtm(:,:,1,band) + alpha_strat(band)*tau_aer_strat(:,:)
173    ENDDO
174
175end subroutine readaerosolstrato_rrtm
Note: See TracBrowser for help on using the repository browser.