source: LMDZ6/branches/Ocean_skin/libf/phylmd/readaerosolstrato.F90 @ 5453

Last change on this file since 5453 was 3605, checked in by lguez, 5 years ago

Merge revisions 3427:3600 of trunk into branch Ocean_skin

  • 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.9 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 print_control_mod, ONLY: prt_level,lunout
18#ifdef CPP_XIOS
19    USE xios
20#endif
21    implicit none
22
23    include "YOMCST.h"
24
25! Variable input
26    logical debut
27
28! Variables locales
29    integer n_lat   ! number of latitudes in the input data
30    integer n_lon   ! number of longitudes in the input data
31    integer n_lev   ! number of levels in the input data
32    integer n_month ! number of months in the input data
33    real, pointer:: latitude(:)
34    real, pointer:: longitude(:)
35    real, pointer:: time(:)
36    real, pointer:: lev(:)
37    integer i, k, band, wave
38    integer, save :: mth_pre=1
39!$OMP THREADPRIVATE(mth_pre)
40
41    real, allocatable, dimension(:,:), save :: tau_aer_strat
42!$OMP THREADPRIVATE(tau_aer_strat)
43
44! Champs reconstitues
45    real, allocatable:: tauaerstrat(:, :, :, :)
46    real, allocatable:: tauaerstrat_mois(:, :, :)
47    real, allocatable:: tauaerstrat_mois_glo(:, :)
48    real, allocatable:: tau_aer_strat_mpi(:, :)
49
50! For NetCDF:
51    integer ncid_in  ! IDs for input files
52    integer varid, ncerr
53
54! Stratospheric aerosols optical properties
55! alpha_strat over the 2 bands is normalised by the 550 nm extinction coefficient
56! alpha_strat_wave is *not* normalised by the 550 nm extinction coefficient
57    real, dimension(nbands) :: alpha_strat, piz_strat, cg_strat
58    data alpha_strat/0.9922547, 0.7114912 /
59    data piz_strat  /0.9999998, 0.99762493/
60    data cg_strat   /0.73107845,0.73229635/
61    real, dimension(nwave_sw) :: alpha_strat_wave
62    data alpha_strat_wave/3.36780953,3.34667683,3.20444202,3.0293026,2.82108808/
63
64    CHARACTER (len = 20)                      :: modname = 'readaerosolstrato'
65    CHARACTER (len = 80)                      :: abort_message
66
67!--------------------------------------------------------
68
69    IF (.not.ALLOCATED(tau_aer_strat)) ALLOCATE(tau_aer_strat(klon,klev))
70
71!--only read file if beginning of run or start of new month
72    IF (debut.OR.mth_cur.NE.mth_pre) THEN
73
74!--only root reads
75    IF (is_mpi_root.AND.is_omp_root) THEN
76
77    IF (nbands.NE.2) THEN
78        abort_message='nbands doit etre egal a 2 dans readaerosolstrat'
79        CALL abort_physic(modname,abort_message,1)
80    ENDIF
81
82    CALL nf95_open("taustrat.nc", nf90_nowrite, ncid_in)
83
84    CALL nf95_inq_varid(ncid_in, "LEV", varid)
85    CALL nf95_gw_var(ncid_in, varid, lev)
86    n_lev = size(lev)
87    IF (n_lev.NE.klev) THEN
88       abort_message='Le nombre de niveaux n est pas egal a klev'
89       CALL abort_physic(modname,abort_message,1)
90    ENDIF
91
92    CALL nf95_inq_varid(ncid_in, "LAT", varid)
93    CALL nf95_gw_var(ncid_in, varid, latitude)
94    n_lat = size(latitude)
95    WRITE(lunout,*) 'LAT aerosol strato=', n_lat, latitude
96    IF (grid_type/=unstructured) THEN
97      IF (n_lat.NE.nbp_lat) THEN
98         abort_message='Le nombre de lat n est pas egal a nbp_lat'
99         CALL abort_physic(modname,abort_message,1)
100      ENDIF
101    ENDIF
102   
103    CALL nf95_inq_varid(ncid_in, "LON", varid)
104    CALL nf95_gw_var(ncid_in, varid, longitude)
105    n_lon = size(longitude)
106    IF (grid_type/=unstructured) THEN
107      WRITE(lunout,*) 'LON aerosol strato=', n_lon, longitude
108      IF (n_lon.NE.nbp_lon) THEN
109         abort_message='Le nombre de lon n est pas egal a nbp_lon'
110         CALL abort_physic(modname,abort_message,1)
111      ENDIF
112    ENDIF
113   
114    CALL nf95_inq_varid(ncid_in, "TIME", varid)
115    CALL nf95_gw_var(ncid_in, varid, time)
116    n_month = size(time)
117    WRITE(lunout,*) 'TIME aerosol strato=', n_month, time
118    IF (n_month.NE.12) THEN
119       abort_message='Le nombre de month n est pas egal a 12'
120       CALL abort_physic(modname,abort_message,1)
121    ENDIF
122
123    IF (.not.ALLOCATED(tauaerstrat))          ALLOCATE(tauaerstrat(n_lon, n_lat, n_lev, n_month))
124    IF (.not.ALLOCATED(tauaerstrat_mois))     ALLOCATE(tauaerstrat_mois(n_lon, n_lat, n_lev))
125    IF (.not.ALLOCATED(tauaerstrat_mois_glo)) ALLOCATE(tauaerstrat_mois_glo(klon_glo, n_lev))
126
127!--reading stratospheric AOD at 550 nm
128    CALL nf95_inq_varid(ncid_in, "TAUSTRAT", varid)
129    ncerr = nf90_get_var(ncid_in, varid, tauaerstrat)
130    WRITE(lunout,*) 'code erreur readaerosolstrato=', ncerr, varid
131
132    CALL nf95_close(ncid_in)
133
134!---select the correct month
135    IF (mth_cur.LT.1.OR.mth_cur.GT.12) THEN
136     WRITE(lunout,*) 'probleme avec le mois dans readaerosolstrat =', mth_cur
137    ENDIF
138    tauaerstrat_mois(:,:,:) = tauaerstrat(:,:,:,mth_cur)
139
140!---reduce to a klon_glo grid
141    CALL grid2dTo1d_glo(tauaerstrat_mois,tauaerstrat_mois_glo)
142
143    ELSE
144      ALLOCATE(tauaerstrat_mois(0,0,0))
145    ENDIF !--is_mpi_root and is_omp_root
146
147!$OMP BARRIER
148
149    IF (grid_type==unstructured) THEN
150#ifdef CPP_XIOS
151      IF (is_omp_master) THEN
152        CALL xios_send_field("taustrat_in",tauaerstrat_mois)
153        ALLOCATE(tau_aer_strat_mpi(klon_mpi, klev))
154        CALL xios_recv_field("taustrat_out",tau_aer_strat_mpi)
155      ELSE
156        ALLOCATE(tau_aer_strat_mpi(0,0))
157      ENDIF
158      CALL scatter_omp(tau_aer_strat_mpi,tau_aer_strat)
159#endif
160    ELSE 
161!--scatter on all proc
162      CALL scatter(tauaerstrat_mois_glo,tau_aer_strat)
163    ENDIF
164
165!--keep memory of previous month
166    mth_pre=mth_cur
167!
168    IF (is_mpi_root.AND.is_omp_root) THEN
169!
170    DEALLOCATE(tauaerstrat)
171    DEALLOCATE(tauaerstrat_mois)
172    DEALLOCATE(tauaerstrat_mois_glo)
173!
174    ENDIF !-is_mpi_root and is_omp_root
175
176!$OMP BARRIER
177
178    ENDIF !--debut ou nouveau mois
179
180!--total vertical aod at the 6 wavelengths
181    DO wave=1, nwave_sw
182    DO k=1, klev
183    tausum_aero(:,wave,id_STRAT_phy)=tausum_aero(:,wave,id_STRAT_phy)+tau_aer_strat(:,k)*alpha_strat_wave(wave)/alpha_strat_wave(2)
184    ENDDO
185    ENDDO
186
187!--weighted average for cg, piz and tau, adding strat aerosols on top of tropospheric ones
188    DO band=1, nbands
189!--anthropogenic aerosols bands 1 and 2
190    cg_aero(:,:,3,band)  = ( cg_aero(:,:,3,band)*piz_aero(:,:,3,band)*tau_aero(:,:,3,band) +          &
191                             cg_strat(band)*piz_strat(band)*alpha_strat(band)*tau_aer_strat(:,:) ) /  &
192                             MAX( piz_aero(:,:,3,band)*tau_aero(:,:,3,band) +                         &
193                                  piz_strat(band)*alpha_strat(band)*tau_aer_strat(:,:), 1.e-15 )
194    piz_aero(:,:,3,band)  = ( piz_aero(:,:,3,band)*tau_aero(:,:,3,band) +                             &
195                              piz_strat(band)*alpha_strat(band)*tau_aer_strat(:,:) ) /                &
196                              MAX( tau_aero(:,:,3,band) + alpha_strat(band)*tau_aer_strat(:,:), 1.e-15 )
197    tau_aero(:,:,3,band)  = tau_aero(:,:,3,band) + alpha_strat(band)*tau_aer_strat(:,:)
198!--natural aerosols bands 1 and 2
199    cg_aero(:,:,2,band)  = ( cg_aero(:,:,2,band)*piz_aero(:,:,2,band)*tau_aero(:,:,2,band) +          &
200                             cg_strat(band)*piz_strat(band)*alpha_strat(band)*tau_aer_strat(:,:) ) /  &
201                             MAX( piz_aero(:,:,2,band)*tau_aero(:,:,2,band) +                         &
202                                  piz_strat(band)*alpha_strat(band)*tau_aer_strat(:,:), 1.e-15 )
203    piz_aero(:,:,2,band)  = ( piz_aero(:,:,2,band)*tau_aero(:,:,2,band) +                             &
204                              piz_strat(band)*alpha_strat(band)*tau_aer_strat(:,:) ) /                &
205                              MAX( tau_aero(:,:,2,band) + alpha_strat(band)*tau_aer_strat(:,:),1.e-15 )
206    tau_aero(:,:,2,band)  = tau_aero(:,:,2,band) + alpha_strat(band)*tau_aer_strat(:,:)
207    ENDDO
208
209end subroutine readaerosolstrato
Note: See TracBrowser for help on using the repository browser.