source: LMDZ5/branches/testing/libf/phylmd/rrtm/readaerosolstrato2_rrtm.F90 @ 2720

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

Merged trunk changes r2664:2719 into testing branch

File size: 12.9 KB
Line 
1!
2! $Id: readaerosolstrato2_rrtm.F90 2526 2016-05-26 22:13:40Z oboucher $
3!
4SUBROUTINE readaerosolstrato2_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, ONLY: nbp_lon, nbp_lat, klon_glo, grid2dTo1d_glo
12    USE mod_phys_lmdz_mpi_data, ONLY :  is_mpi_root
13    USE mod_phys_lmdz_omp_data, ONLY :  is_omp_root
14    USE mod_phys_lmdz_para
15    USE phys_state_var_mod
16    USE phys_local_var_mod
17    USE aero_mod
18    USE dimphy
19    USE YOERAD, ONLY : NLW
20    USE YOMCST
21
22    IMPLICIT NONE
23
24    INCLUDE "clesphys.h"
25
26    CHARACTER (len = 80) :: abort_message
27    CHARACTER (LEN=20) :: modname = 'readaerosolstrato2'
28
29! Variable input
30    LOGICAL, INTENT(IN) ::  debut
31
32! Variables locales
33    INTEGER n_lat   ! number of latitudes in the input data
34    INTEGER n_lon   ! number of longitudes
35    INTEGER n_lev   ! number of levels in the input data
36    INTEGER n_month ! number of months in the input data
37    INTEGER n_wav   ! number of wavelengths in the input data
38    REAL, POINTER:: latitude(:)
39    REAL, POINTER:: time(:)
40    REAL, POINTER:: lev(:)
41    REAL, POINTER:: wav(:)
42    INTEGER i,k,wave,band
43    INTEGER, SAVE :: mth_pre
44
45    REAL, ALLOCATABLE, DIMENSION(:,:,:), SAVE :: tau_aer_strat
46    REAL, ALLOCATABLE, DIMENSION(:,:,:), SAVE :: piz_aer_strat
47    REAL, ALLOCATABLE, DIMENSION(:,:,:), SAVE :: cg_aer_strat
48    REAL, ALLOCATABLE, DIMENSION(:,:,:), SAVE :: taulw_aer_strat
49!$OMP THREADPRIVATE(tau_aer_strat,piz_aer_strat,cg_aer_strat,taulw_aer_strat)
50
51! Champs reconstitues
52    REAL, ALLOCATABLE:: tauaerstrat(:, :, :, :)
53    REAL, ALLOCATABLE:: pizaerstrat(:, :, :, :)
54    REAL, ALLOCATABLE:: cgaerstrat(:, :, :, :)
55    REAL, ALLOCATABLE:: taulwaerstrat(:, :, :, :)
56
57    REAL, ALLOCATABLE:: tauaerstrat_mois(:, :, :, :)
58    REAL, ALLOCATABLE:: pizaerstrat_mois(:, :, :, :)
59    REAL, ALLOCATABLE:: cgaerstrat_mois(:, :, :, :)
60    REAL, ALLOCATABLE:: taulwaerstrat_mois(:, :, :, :)
61
62    REAL, ALLOCATABLE:: tauaerstrat_mois_glo(:, :, :)
63    REAL, ALLOCATABLE:: pizaerstrat_mois_glo(:, :, :)
64    REAL, ALLOCATABLE:: cgaerstrat_mois_glo(:, :, :)
65    REAL, ALLOCATABLE:: taulwaerstrat_mois_glo(:, :, :)
66
67! For NetCDF:
68    INTEGER ncid_in  ! IDs for input files
69    INTEGER varid, ncerr
70
71!--------------------------------------------------------
72
73    IF (.not.ALLOCATED(tau_aer_strat)) ALLOCATE(tau_aer_strat(klon,klev,NSW))
74    IF (.not.ALLOCATED(piz_aer_strat)) ALLOCATE(piz_aer_strat(klon,klev,NSW))
75    IF (.not.ALLOCATED(cg_aer_strat))  ALLOCATE(cg_aer_strat(klon,klev,NSW))
76
77    IF (.not.ALLOCATED(taulw_aer_strat)) ALLOCATE(taulw_aer_strat(klon,klev,NLW))
78
79!--we only read monthly strat aerosol data
80    IF (debut.OR.mth_cur.NE.mth_pre) THEN
81
82!--only root reads the data
83    IF (is_mpi_root.AND.is_omp_root) THEN
84
85!--check mth_cur
86    IF (mth_cur.LT.1.OR.mth_cur.GT.12) THEN
87      print *,'probleme avec le mois dans readaerosolstrat =', mth_cur
88    ENDIF
89
90!--initialize n_lon as input data is 2D (lat-alt) only
91    n_lon = nbp_lon
92
93!--Starts with SW optical properties
94
95    CALL nf95_open("tauswstrat.2D.nc", nf90_nowrite, ncid_in)
96
97    CALL nf95_inq_varid(ncid_in, "LEV", varid)
98    CALL nf95_gw_var(ncid_in, varid, lev)
99    n_lev = size(lev)
100    IF (n_lev.NE.klev) THEN
101       abort_message='Le nombre de niveaux n est pas egal a klev'
102       CALL abort_physic(modname,abort_message,1)
103    ENDIF
104
105    CALL nf95_inq_varid(ncid_in, "LAT", varid)
106    CALL nf95_gw_var(ncid_in, varid, latitude)
107    n_lat = size(latitude)
108    IF (n_lat.NE.nbp_lat) THEN
109       print *, 'latitude=', n_lat, nbp_lat
110       abort_message='Le nombre de lat n est pas egal a nbp_lat'
111       CALL abort_physic(modname,abort_message,1)
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    IF (n_month.NE.12) THEN
118       abort_message='Le nombre de month n est pas egal a 12'
119       CALL abort_physic(modname,abort_message,1)
120    ENDIF
121
122    CALL nf95_inq_varid(ncid_in, "WAV", varid)
123    CALL nf95_gw_var(ncid_in, varid, wav)
124    n_wav = size(wav)
125    print *, 'WAV aerosol strato=', n_wav, wav
126    IF (n_wav.NE.NSW) THEN
127       abort_message='Le nombre de wav n est pas egal a NSW'
128       CALL abort_physic(modname,abort_message,1)
129    ENDIF
130
131    ALLOCATE(tauaerstrat(n_lat, n_lev, n_wav, n_month))
132    ALLOCATE(pizaerstrat(n_lat, n_lev, n_wav, n_month))
133    ALLOCATE(cgaerstrat(n_lat, n_lev, n_wav, n_month))
134
135    ALLOCATE(tauaerstrat_mois(n_lon, n_lat, n_lev, n_wav))
136    ALLOCATE(pizaerstrat_mois(n_lon, n_lat, n_lev, n_wav))
137    ALLOCATE(cgaerstrat_mois(n_lon, n_lat, n_lev, n_wav))
138
139    ALLOCATE(tauaerstrat_mois_glo(klon_glo, n_lev, n_wav))
140    ALLOCATE(pizaerstrat_mois_glo(klon_glo, n_lev, n_wav))
141    ALLOCATE(cgaerstrat_mois_glo(klon_glo, n_lev, n_wav))
142
143!--reading stratospheric aerosol tau per layer
144    CALL nf95_inq_varid(ncid_in, "TAU_SUN", varid)
145    ncerr = nf90_get_var(ncid_in, varid, tauaerstrat)
146    print *,'code erreur readaerosolstrato=', ncerr, varid
147
148!--reading stratospheric aerosol omega per layer
149    CALL nf95_inq_varid(ncid_in, "OME_SUN", varid)
150    ncerr = nf90_get_var(ncid_in, varid, pizaerstrat)
151    print *,'code erreur readaerosolstrato=', ncerr, varid
152
153!--reading stratospheric aerosol g per layer
154    CALL nf95_inq_varid(ncid_in, "GGG_SUN", varid)
155    ncerr = nf90_get_var(ncid_in, varid, cgaerstrat)
156    print *,'code erreur readaerosolstrato sw=', ncerr, varid
157
158    CALL nf95_close(ncid_in)
159
160!--select the correct month
161!--and copy into 1st longitude
162    tauaerstrat_mois(1,:,:,:) = tauaerstrat(:,:,:,mth_cur)
163    pizaerstrat_mois(1,:,:,:) = pizaerstrat(:,:,:,mth_cur)
164    cgaerstrat_mois(1,:,:,:)  = cgaerstrat(:,:,:,mth_cur)
165
166!--copy longitudes
167    DO i=2, n_lon
168     tauaerstrat_mois(i,:,:,:) = tauaerstrat_mois(1,:,:,:)
169     pizaerstrat_mois(i,:,:,:) = pizaerstrat_mois(1,:,:,:)
170      cgaerstrat_mois(i,:,:,:)  = cgaerstrat_mois(1,:,:,:)
171    ENDDO
172
173!---reduce to a klon_glo grid
174    DO band=1, NSW
175      CALL grid2dTo1d_glo(tauaerstrat_mois(:,:,:,band),tauaerstrat_mois_glo(:,:,band))
176      CALL grid2dTo1d_glo(pizaerstrat_mois(:,:,:,band),pizaerstrat_mois_glo(:,:,band))
177      CALL grid2dTo1d_glo(cgaerstrat_mois(:,:,:,band),cgaerstrat_mois_glo(:,:,band))
178    ENDDO
179
180!--Now LW optical properties
181!
182    CALL nf95_open("taulwstrat.2D.nc", nf90_nowrite, ncid_in)
183
184    CALL nf95_inq_varid(ncid_in, "LEV", varid)
185    CALL nf95_gw_var(ncid_in, varid, lev)
186    n_lev = size(lev)
187    IF (n_lev.NE.klev) THEN
188       abort_message='Le nombre de niveaux n est pas egal a klev'
189       CALL abort_physic(modname,abort_message,1)
190    ENDIF
191
192    CALL nf95_inq_varid(ncid_in, "LAT", varid)
193    CALL nf95_gw_var(ncid_in, varid, latitude)
194    n_lat = size(latitude)
195    IF (n_lat.NE.nbp_lat) THEN
196       abort_message='Le nombre de lat n est pas egal a nbp_lat'
197       CALL abort_physic(modname,abort_message,1)
198    ENDIF
199
200    CALL nf95_inq_varid(ncid_in, "TIME", varid)
201    CALL nf95_gw_var(ncid_in, varid, time)
202    n_month = size(time)
203    IF (n_month.NE.12) THEN
204       abort_message='Le nombre de month n est pas egal a 12'
205       CALL abort_physic(modname,abort_message,1)
206    ENDIF
207
208    CALL nf95_inq_varid(ncid_in, "WAV", varid)
209    CALL nf95_gw_var(ncid_in, varid, wav)
210    n_wav = size(wav)
211    print *, 'WAV aerosol strato=', n_wav, wav
212    IF (n_wav.NE.NLW) THEN
213       abort_message='Le nombre de wav n est pas egal a NLW'
214       CALL abort_physic(modname,abort_message,1)
215    ENDIF
216
217    ALLOCATE(taulwaerstrat(n_lat, n_lev, n_wav, n_month))
218    ALLOCATE(taulwaerstrat_mois(n_lon, n_lat, n_lev, n_wav))
219    ALLOCATE(taulwaerstrat_mois_glo(klon_glo, n_lev, n_wav))
220
221!--reading stratospheric aerosol lw tau per layer
222    CALL nf95_inq_varid(ncid_in, "TAU_EAR", varid)
223    ncerr = nf90_get_var(ncid_in, varid, taulwaerstrat)
224    print *,'code erreur readaerosolstrato lw=', ncerr, varid
225
226    CALL nf95_close(ncid_in)
227
228!--select the correct month
229!--and copy into 1st longitude
230    taulwaerstrat_mois(1,:,:,:) = taulwaerstrat(:,:,:,mth_cur)
231!--copy longitudes
232    DO i=2, n_lon
233      taulwaerstrat_mois(i,:,:,:) = taulwaerstrat_mois(1,:,:,:)
234    ENDDO
235
236!---reduce to a klon_glo grid
237    DO band=1, NLW
238      CALL grid2dTo1d_glo(taulwaerstrat_mois(:,:,:,band),taulwaerstrat_mois_glo(:,:,band))
239    ENDDO
240
241    ENDIF !--is_mpi_root and is_omp_root
242
243!$OMP BARRIER
244
245!--keep memory of previous month
246    mth_pre=mth_cur
247
248!--scatter on all proc
249    CALL scatter(tauaerstrat_mois_glo,tau_aer_strat)
250    CALL scatter(pizaerstrat_mois_glo,piz_aer_strat)
251    CALL scatter(cgaerstrat_mois_glo,cg_aer_strat)
252    CALL scatter(taulwaerstrat_mois_glo,taulw_aer_strat)
253
254    IF (is_mpi_root.AND.is_omp_root) THEN
255!
256    DEALLOCATE(tauaerstrat, pizaerstrat, cgaerstrat)
257    DEALLOCATE(tauaerstrat_mois, pizaerstrat_mois, cgaerstrat_mois)
258    DEALLOCATE(tauaerstrat_mois_glo,pizaerstrat_mois_glo,cgaerstrat_mois_glo)
259
260    DEALLOCATE(taulwaerstrat,taulwaerstrat_mois,taulwaerstrat_mois_glo)
261!
262    ENDIF !--is_mpi_root and is_omp_root
263
264!$OMP BARRIER
265
266    ENDIF !--debut ou nouveau mois
267
268!--total vertical aod at the 5 SW wavelengths
269!--for now use band 3 AOD into all 5 wavelengths
270!--it is only a reasonable approximation for 550 nm (wave=2)
271    band=3
272    DO i=1, klon
273    DO k=1, klev
274      IF (stratomask(i,k).GT.0.999999) THEN
275        DO wave=1, nwave_sw
276          tausum_aero(i,wave,id_STRAT_phy)=tausum_aero(i,wave,id_STRAT_phy)+tau_aer_strat(i,k,band)
277        ENDDO
278      ENDIF
279    ENDDO
280    ENDDO
281
282!--weighted average for cg, piz and tau, adding strat aerosols on top of tropospheric ones
283    DO band=1, NSW
284    WHERE (stratomask.GT.0.999999)
285!--anthropogenic aerosols bands 1 to NSW
286    cg_aero_sw_rrtm(:,:,2,band)  = ( cg_aero_sw_rrtm(:,:,2,band)*piz_aero_sw_rrtm(:,:,2,band)*tau_aero_sw_rrtm(:,:,2,band) + &
287                                     cg_aer_strat(:,:,band)*piz_aer_strat(:,:,band)*tau_aer_strat(:,:,band) ) /              &
288                                MAX( piz_aero_sw_rrtm(:,:,2,band)*tau_aero_sw_rrtm(:,:,2,band) +                             &
289                                     piz_aer_strat(:,:,band)*tau_aer_strat(:,:,band), 1.e-15 )
290    piz_aero_sw_rrtm(:,:,2,band) = ( piz_aero_sw_rrtm(:,:,2,band)*tau_aero_sw_rrtm(:,:,2,band) +                             &
291                                     piz_aer_strat(:,:,band)*tau_aer_strat(:,:,band) ) /                                     &
292                                MAX( tau_aero_sw_rrtm(:,:,2,band) + tau_aer_strat(:,:,band), 1.e-15 )
293    tau_aero_sw_rrtm(:,:,2,band)  = tau_aero_sw_rrtm(:,:,2,band) + tau_aer_strat(:,:,band)
294!--natural aerosols bands 1 to NSW
295    cg_aero_sw_rrtm(:,:,1,band)  = ( cg_aero_sw_rrtm(:,:,1,band)*piz_aero_sw_rrtm(:,:,1,band)*tau_aero_sw_rrtm(:,:,1,band) + &
296                                     cg_aer_strat(:,:,band)*piz_aer_strat(:,:,band)*tau_aer_strat(:,:,band) ) /              &
297                                MAX( piz_aero_sw_rrtm(:,:,1,band)*tau_aero_sw_rrtm(:,:,1,band) +                             &
298                                     piz_aer_strat(:,:,band)*tau_aer_strat(:,:,band), 1.e-15 )
299    piz_aero_sw_rrtm(:,:,1,band) = ( piz_aero_sw_rrtm(:,:,1,band)*tau_aero_sw_rrtm(:,:,1,band) +                             &
300                                     piz_aer_strat(:,:,band)*tau_aer_strat(:,:,band) ) /                                     &
301                                MAX( tau_aero_sw_rrtm(:,:,1,band) + tau_aer_strat(:,:,band), 1.e-15 )
302    tau_aero_sw_rrtm(:,:,1,band)  = tau_aero_sw_rrtm(:,:,1,band) + tau_aer_strat(:,:,band)
303!--no stratospheric aerosol in index 1 for these tests
304!    cg_aero_sw_rrtm(:,:,1,band)  =  cg_aero_sw_rrtm(:,:,1,band)
305!    piz_aero_sw_rrtm(:,:,1,band)  = piz_aero_sw_rrtm(:,:,1,band)
306!    tau_aero_sw_rrtm(:,:,1,band)  = tau_aero_sw_rrtm(:,:,1,band)
307    ENDWHERE
308    ENDDO
309
310!--total vertical aod at 10 um
311!--this is approximated from band 7 of RRTM
312    band=7
313    DO i=1, klon
314    DO k=1, klev
315      IF (stratomask(i,k).GT.0.999999) THEN
316        DO wave=1, nwave_lw
317          tausum_aero(i,nwave_sw+wave,id_STRAT_phy)=tausum_aero(i,nwave_sw+wave,id_STRAT_phy)+taulw_aer_strat(i,k,band)
318        ENDDO
319      ENDIF
320    ENDDO
321    ENDDO
322
323    DO band=1, NLW
324    WHERE (stratomask.GT.0.999999)
325    tau_aero_lw_rrtm(:,:,2,band)  = tau_aero_lw_rrtm(:,:,2,band) + taulw_aer_strat(:,:,band)
326    tau_aero_lw_rrtm(:,:,1,band)  = tau_aero_lw_rrtm(:,:,1,band) + taulw_aer_strat(:,:,band)
327!--no stratospheric aerosols in index 1 for these tests
328!    tau_aero_lw_rrtm(:,:,1,band)  = tau_aero_lw_rrtm(:,:,1,band)
329    ENDWHERE
330    ENDDO
331
332!--default SSA value if there is no aerosol
333!--to avoid 0 values that seems to cause some problem to RRTM
334    WHERE (tau_aero_sw_rrtm.LT.1.e-14)
335      piz_aero_sw_rrtm = 1.0
336    ENDWHERE
337
338!--in principle this should not be necessary
339!--as these variables have min values already but just in case
340!--put 1e-15 min value to both SW and LW AOD
341    tau_aero_sw_rrtm = MAX(tau_aero_sw_rrtm,1.e-15)
342    tau_aero_lw_rrtm = MAX(tau_aero_lw_rrtm,1.e-15)
343
344END SUBROUTINE readaerosolstrato2_rrtm
Note: See TracBrowser for help on using the repository browser.