source: LMDZ5/trunk/libf/phylmd/rrtm/readaerosolstrato_rrtm.F90 @ 2058

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

Modifications sur l'utilisation du module/fichier YOMCST


Modification on the use of the file/module YOMCST
MPL

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