source: LMDZ6/branches/DYNAMICO-conv/libf/phylmd/readaerosolstrato.F90 @ 3863

Last change on this file since 3863 was 3413, checked in by Laurent Fairhead, 6 years ago

Inclusion of Yann's latest (summer/fall 2018) modifications for
convergence of DYNAMICO/LMDZ physics
YM/LF

  • 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: 7.4 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, ONLY: nbp_lon, nbp_lat, klon_glo, &
9                                 grid2dto1d_glo, grid_type, unstructured
10    USE mod_phys_lmdz_mpi_data, ONLY :  is_mpi_root
11    USE mod_phys_lmdz_omp_data, ONLY :  is_omp_root
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    USE xios
18
19    implicit none
20
21    include "YOMCST.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 i, k, band, wave
36    integer, save :: mth_pre=1
37!$OMP THREADPRIVATE(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:: tau_aer_strat_mpi(:, :)
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) :: alpha_strat, piz_strat, cg_strat
56    data alpha_strat/0.9922547, 0.7114912 /
57    data piz_strat  /0.9999998, 0.99762493/
58    data cg_strat   /0.73107845,0.73229635/
59    real, dimension(nwave_sw) :: 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!--only read file if beginning of run or start of new month
67    IF (debut.OR.mth_cur.NE.mth_pre) THEN
68
69!--only root reads
70    IF (is_mpi_root.AND.is_omp_root) THEN
71
72    IF (nbands.NE.2) THEN
73        print *,'nbands doit etre egal a 2 dans readaerosolstrat'
74        STOP
75    ENDIF
76
77    CALL nf95_open("taustrat.nc", nf90_nowrite, ncid_in)
78
79    CALL nf95_inq_varid(ncid_in, "LEV", varid)
80    CALL nf95_gw_var(ncid_in, varid, lev)
81    n_lev = size(lev)
82    IF (n_lev.NE.klev) THEN
83       print *,'Le nombre de niveaux n est pas egal a klev'
84       STOP
85    ENDIF
86
87    CALL nf95_inq_varid(ncid_in, "LAT", varid)
88    CALL nf95_gw_var(ncid_in, varid, latitude)
89    n_lat = size(latitude)
90    print *, 'LAT aerosol strato=', n_lat, latitude
91    IF (grid_type/=unstructured) THEN
92      IF (n_lat.NE.nbp_lat) THEN
93         print *,'Le nombre de lat n est pas egal a nbp_lat'
94         STOP
95      ENDIF
96    ENDIF
97   
98    CALL nf95_inq_varid(ncid_in, "LON", varid)
99    CALL nf95_gw_var(ncid_in, varid, longitude)
100    n_lon = size(longitude)
101    IF (grid_type/=unstructured) THEN
102      print *, 'LON aerosol strato=', n_lon, longitude
103      IF (n_lon.NE.nbp_lon) THEN
104         print *,'Le nombre de lon n est pas egal a nbp_lon'
105         STOP
106      ENDIF
107    ENDIF
108   
109    CALL nf95_inq_varid(ncid_in, "TIME", varid)
110    CALL nf95_gw_var(ncid_in, varid, time)
111    n_month = size(time)
112    print *, 'TIME aerosol strato=', n_month, time
113    IF (n_month.NE.12) THEN
114       print *,'Le nombre de month n est pas egal a 12'
115       STOP
116    ENDIF
117
118    IF (.not.ALLOCATED(tauaerstrat))          ALLOCATE(tauaerstrat(n_lon, n_lat, n_lev, n_month))
119    IF (.not.ALLOCATED(tauaerstrat_mois))     ALLOCATE(tauaerstrat_mois(n_lon, n_lat, n_lev))
120    IF (.not.ALLOCATED(tauaerstrat_mois_glo)) ALLOCATE(tauaerstrat_mois_glo(klon_glo, n_lev))
121
122!--reading stratospheric AOD at 550 nm
123    CALL nf95_inq_varid(ncid_in, "TAUSTRAT", varid)
124    ncerr = nf90_get_var(ncid_in, varid, tauaerstrat)
125    print *,'code erreur readaerosolstrato=', ncerr, varid
126
127    CALL nf95_close(ncid_in)
128
129!---select the correct month
130    IF (mth_cur.LT.1.OR.mth_cur.GT.12) THEN
131      print *,'probleme avec le mois dans readaerosolstrat =', mth_cur
132    ENDIF
133    tauaerstrat_mois(:,:,:) = tauaerstrat(:,:,:,mth_cur)
134
135!---reduce to a klon_glo grid
136    CALL grid2dTo1d_glo(tauaerstrat_mois,tauaerstrat_mois_glo)
137
138    ELSE
139      ALLOCATE(tauaerstrat_mois(0,0,0))
140    ENDIF !--is_mpi_root and is_omp_root
141
142!$OMP BARRIER
143
144    IF (grid_type==unstructured) THEN
145      IF (is_omp_master) THEN
146        CALL xios_send_field("taustrat_in",tauaerstrat_mois)
147        ALLOCATE(tau_aer_strat_mpi(klon_mpi, klev))
148        CALL xios_recv_field("taustrat_out",tau_aer_strat_mpi)
149      ELSE
150        ALLOCATE(tau_aer_strat_mpi(0,0))
151      ENDIF
152      CALL scatter_omp(tau_aer_strat_mpi,tau_aer_strat)
153    ELSE 
154!--scatter on all proc
155      CALL scatter(tauaerstrat_mois_glo,tau_aer_strat)
156    ENDIF
157
158!--keep memory of previous month
159    mth_pre=mth_cur
160!
161    IF (is_mpi_root.AND.is_omp_root) THEN
162!
163    DEALLOCATE(tauaerstrat)
164    DEALLOCATE(tauaerstrat_mois)
165    DEALLOCATE(tauaerstrat_mois_glo)
166!
167    ENDIF !-is_mpi_root and is_omp_root
168
169!$OMP BARRIER
170
171    ENDIF !--debut ou nouveau mois
172
173!--total vertical aod at the 6 wavelengths
174    DO wave=1, nwave_sw
175    DO k=1, klev
176    tausum_aero(:,wave,id_STRAT_phy)=tausum_aero(:,wave,id_STRAT_phy)+tau_aer_strat(:,k)*alpha_strat_wave(wave)/alpha_strat_wave(2)
177    ENDDO
178    ENDDO
179
180!--weighted average for cg, piz and tau, adding strat aerosols on top of tropospheric ones
181    DO band=1, nbands
182!--anthropogenic aerosols bands 1 and 2
183    cg_aero(:,:,3,band)  = ( cg_aero(:,:,3,band)*piz_aero(:,:,3,band)*tau_aero(:,:,3,band) +          &
184                             cg_strat(band)*piz_strat(band)*alpha_strat(band)*tau_aer_strat(:,:) ) /  &
185                             MAX( piz_aero(:,:,3,band)*tau_aero(:,:,3,band) +                         &
186                                  piz_strat(band)*alpha_strat(band)*tau_aer_strat(:,:), 1.e-15 )
187    piz_aero(:,:,3,band)  = ( piz_aero(:,:,3,band)*tau_aero(:,:,3,band) +                             &
188                              piz_strat(band)*alpha_strat(band)*tau_aer_strat(:,:) ) /                &
189                              MAX( tau_aero(:,:,3,band) + alpha_strat(band)*tau_aer_strat(:,:), 1.e-15 )
190    tau_aero(:,:,3,band)  = tau_aero(:,:,3,band) + alpha_strat(band)*tau_aer_strat(:,:)
191!--natural aerosols bands 1 and 2
192    cg_aero(:,:,2,band)  = ( cg_aero(:,:,2,band)*piz_aero(:,:,2,band)*tau_aero(:,:,2,band) +          &
193                             cg_strat(band)*piz_strat(band)*alpha_strat(band)*tau_aer_strat(:,:) ) /  &
194                             MAX( piz_aero(:,:,2,band)*tau_aero(:,:,2,band) +                         &
195                                  piz_strat(band)*alpha_strat(band)*tau_aer_strat(:,:), 1.e-15 )
196    piz_aero(:,:,2,band)  = ( piz_aero(:,:,2,band)*tau_aero(:,:,2,band) +                             &
197                              piz_strat(band)*alpha_strat(band)*tau_aer_strat(:,:) ) /                &
198                              MAX( tau_aero(:,:,2,band) + alpha_strat(band)*tau_aer_strat(:,:),1.e-15 )
199    tau_aero(:,:,2,band)  = tau_aero(:,:,2,band) + alpha_strat(band)*tau_aer_strat(:,:)
200    ENDDO
201
202end subroutine readaerosolstrato
Note: See TracBrowser for help on using the repository browser.