source: LMDZ5/branches/testing/libf/phylmd/readaerosolstrato.F90 @ 5423

Last change on this file since 5423 was 2787, checked in by Laurent Fairhead, 8 years ago

Merged trunk changes r2727:2785 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
File size: 6.7 KB
RevLine 
[1764]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
[2408]8    USE mod_grid_phy_lmdz, ONLY: nbp_lon, nbp_lat, klon_glo, &
9                                 grid2dto1d_glo
[1764]10    USE mod_phys_lmdz_mpi_data, ONLY :  is_mpi_root
[2542]11    USE mod_phys_lmdz_omp_data, ONLY :  is_omp_root
[1764]12    USE mod_phys_lmdz_para
13    USE phys_state_var_mod
14    USE phys_local_var_mod
15    USE aero_mod
16    USE dimphy
17
18    implicit none
19
20    include "YOMCST.h"
21
22! Variable input
23    logical debut
24
25! Variables locales
26    integer n_lat   ! number of latitudes in the input data
27    integer n_lon   ! number of longitudes in the input data
28    integer n_lev   ! number of levels in the input data
29    integer n_month ! number of months in the input data
30    real, pointer:: latitude(:)
31    real, pointer:: longitude(:)
32    real, pointer:: time(:)
33    real, pointer:: lev(:)
34    integer i, k, band, wave
[2787]35    integer, save :: mth_pre=1
[1864]36!$OMP THREADPRIVATE(mth_pre)
[1764]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
46! For NetCDF:
47    integer ncid_in  ! IDs for input files
48    integer varid, ncerr
49
50! Stratospheric aerosols optical properties
51! alpha_strat over the 2 bands is normalised by the 550 nm extinction coefficient
52! alpha_strat_wave is *not* normalised by the 550 nm extinction coefficient
53    real, dimension(nbands) :: alpha_strat, piz_strat, cg_strat
54    data alpha_strat/0.9922547, 0.7114912 /
55    data piz_strat  /0.9999998, 0.99762493/
56    data cg_strat   /0.73107845,0.73229635/
[2594]57    real, dimension(nwave_sw) :: alpha_strat_wave
[1764]58    data alpha_strat_wave/3.36780953,3.34667683,3.20444202,3.0293026,2.82108808/
59
60!--------------------------------------------------------
61
62    IF (.not.ALLOCATED(tau_aer_strat)) ALLOCATE(tau_aer_strat(klon,klev))
63
[1864]64!--only read file if beginning of run or start of new month
65    IF (debut.OR.mth_cur.NE.mth_pre) THEN
66
[2542]67!--only root reads
68    IF (is_mpi_root.AND.is_omp_root) THEN
[1764]69
70    IF (nbands.NE.2) THEN
71        print *,'nbands doit etre egal a 2 dans readaerosolstrat'
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
[2408]89    IF (n_lat.NE.nbp_lat) THEN
90       print *,'Le nombre de lat n est pas egal a nbp_lat'
[1764]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
[2408]98    IF (n_lon.NE.nbp_lon) THEN
99       print *,'Le nombre de lon n est pas egal a nbp_lon'
[1764]100       STOP
101    ENDIF
102
[1771]103    CALL nf95_inq_varid(ncid_in, "TIME", varid)
[1764]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
[1864]112    IF (.not.ALLOCATED(tauaerstrat))          ALLOCATE(tauaerstrat(n_lon, n_lat, n_lev, n_month))
113    IF (.not.ALLOCATED(tauaerstrat_mois))     ALLOCATE(tauaerstrat_mois(n_lon, n_lat, n_lev))
114    IF (.not.ALLOCATED(tauaerstrat_mois_glo)) ALLOCATE(tauaerstrat_mois_glo(klon_glo, n_lev))
[1764]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
[2542]132    ENDIF !--is_mpi_root and is_omp_root
[1864]133
[2542]134!$OMP BARRIER
135
[1764]136!--scatter on all proc
137    CALL scatter(tauaerstrat_mois_glo,tau_aer_strat)
138
[1864]139!--keep memory of previous month
[1764]140    mth_pre=mth_cur
[2160]141!
[2542]142    IF (is_mpi_root.AND.is_omp_root) THEN
[2160]143!
144    DEALLOCATE(tauaerstrat)
145    DEALLOCATE(tauaerstrat_mois)
146    DEALLOCATE(tauaerstrat_mois_glo)
147!
[2542]148    ENDIF !-is_mpi_root and is_omp_root
[1764]149
[2542]150!$OMP BARRIER
151
[1764]152    ENDIF !--debut ou nouveau mois
153
[2160]154!--total vertical aod at the 6 wavelengths
[2594]155    DO wave=1, nwave_sw
[1764]156    DO k=1, klev
[2160]157    tausum_aero(:,wave,id_STRAT_phy)=tausum_aero(:,wave,id_STRAT_phy)+tau_aer_strat(:,k)*alpha_strat_wave(wave)/alpha_strat_wave(2)
[1764]158    ENDDO
159    ENDDO
160
161!--weighted average for cg, piz and tau, adding strat aerosols on top of tropospheric ones
162    DO band=1, nbands
163!--anthropogenic aerosols bands 1 and 2
164    cg_aero(:,:,3,band)  = ( cg_aero(:,:,3,band)*piz_aero(:,:,3,band)*tau_aero(:,:,3,band) +          &
165                             cg_strat(band)*piz_strat(band)*alpha_strat(band)*tau_aer_strat(:,:) ) /  &
166                             MAX( piz_aero(:,:,3,band)*tau_aero(:,:,3,band) +                         &
167                                  piz_strat(band)*alpha_strat(band)*tau_aer_strat(:,:), 1.e-15 )
168    piz_aero(:,:,3,band)  = ( piz_aero(:,:,3,band)*tau_aero(:,:,3,band) +                             &
169                              piz_strat(band)*alpha_strat(band)*tau_aer_strat(:,:) ) /                &
170                              MAX( tau_aero(:,:,3,band) + alpha_strat(band)*tau_aer_strat(:,:), 1.e-15 )
171    tau_aero(:,:,3,band)  = tau_aero(:,:,3,band) + alpha_strat(band)*tau_aer_strat(:,:)
172!--natural aerosols bands 1 and 2
173    cg_aero(:,:,2,band)  = ( cg_aero(:,:,2,band)*piz_aero(:,:,2,band)*tau_aero(:,:,2,band) +          &
174                             cg_strat(band)*piz_strat(band)*alpha_strat(band)*tau_aer_strat(:,:) ) /  &
175                             MAX( piz_aero(:,:,2,band)*tau_aero(:,:,2,band) +                         &
176                                  piz_strat(band)*alpha_strat(band)*tau_aer_strat(:,:), 1.e-15 )
177    piz_aero(:,:,2,band)  = ( piz_aero(:,:,2,band)*tau_aero(:,:,2,band) +                             &
178                              piz_strat(band)*alpha_strat(band)*tau_aer_strat(:,:) ) /                &
179                              MAX( tau_aero(:,:,2,band) + alpha_strat(band)*tau_aer_strat(:,:),1.e-15 )
180    tau_aero(:,:,2,band)  = tau_aero(:,:,2,band) + alpha_strat(band)*tau_aer_strat(:,:)
181    ENDDO
182
183end subroutine readaerosolstrato
Note: See TracBrowser for help on using the repository browser.