source: LMDZ6/branches/LMDZ-INCA-Dyn/libf/phylmd/Dust/precuremission.F @ 5420

Last change on this file since 5420 was 2630, checked in by fhourdin, 8 years ago

Importation du modèle d'aérosols de Boucher, Escribano et al.

File size: 11.8 KB
Line 
1C Subroutine that calculates the emission of aerosols precursors
2      SUBROUTINE precuremission(ftsol,u10m_ec,v10m_ec,
3     .                          pctsrf,u_seri,v_seri,paprs,pplay,cdragh,
4     .                          cdragm,t_seri,q_seri,tsol,fracso2emis,
5     .                          frach2sofso2,bateau,zdz,zalt,
6     .                          kminbc,kmaxbc,pdtphys,scale_param_bb,
7     .                          scale_param_ind,iregion_ind,iregion_bb,
8     .                          nbreg_ind, nbreg_bb,
9     .                          lmt_so2ff_l,lmt_so2ff_h,lmt_so2nff,
10     .                          lmt_so2ba,lmt_so2bb_l,lmt_so2bb_h,
11     .                          lmt_so2volc_cont,lmt_altvolc_cont,
12     .                          lmt_so2volc_expl,lmt_altvolc_expl,
13     .                          lmt_dmsbio,lmt_h2sbio, lmt_dmsconc,
14     .                          lmt_dms,id_prec,id_fine,
15     .                                 flux_sparam_ind,flux_sparam_bb,
16     .                                 source_tr,flux_tr,tr_seri)
17
18      USE dimphy
19      USE indice_sol_mod
20      USE infotrac
21!      USE phytracr_spl_mod, ONLY : nbreg_dust, nbreg_ind, nbreg_bb
22      IMPLICIT NONE
23
24#include "dimensions.h"
25#include "chem.h"
26#include "chem_spla.h"
27#include "YOMCST.h"
28#include "paramet.h"
29
30c============================= INPUT ===================================
31      INTEGER kminbc, kmaxbc
32      REAL ftsol(klon,nbsrf)  ! temperature du sol par type
33      REAL tsol(klon)         ! temperature du sol moyenne
34      REAL t_seri(klon,klev)  ! temperature
35      REAL u_seri(klon,klev)  ! vent
36      REAL v_seri(klon,klev)  ! vent
37      REAL q_seri(klon,klev)  ! vapeur d eau kg/kg
38      REAL u10m_ec(klon), v10m_ec(klon)  ! vent a 10 metres
39      REAL pctsrf(klon,nbsrf)
40      REAL pdtphys  ! pas d'integration pour la physique (seconde)
41      REAL paprs(klon,klev+1)  ! pression pour chaque inter-couche (en Pa)
42      REAL pplay(klon,klev)  ! pression pour le mileu de chaque couche (en Pa)
43      REAL cdragh(klon), cdragm(klon)     
44      REAL fracso2emis        !--fraction so2 emis en so2
45      REAL frach2sofso2       !--fraction h2s from so2
46      REAL zdz(klon,klev)
47      LOGICAL edgar, bateau
48      INTEGER id_prec,id_fine
49c
50c------------------------- Scaling Parameters --------------------------
51c
52      INTEGER nbreg_ind, nbreg_bb
53      INTEGER iregion_ind(klon)  !Defines regions for SO2, BC & OM
54      INTEGER iregion_bb(klon)  !Defines regions for SO2, BC & OM
55      REAL scale_param_bb(nbreg_bb)  !Scaling parameter for biomas burning
56      REAL scale_param_ind(nbreg_ind) !Scaling parameter for industrial emissions
57c
58c============================= OUTPUT ==================================
59c
60      REAL source_tr(klon,nbtr)
61      REAL flux_tr(klon,nbtr)
62      REAL tr_seri(klon,klev,nbtr) ! traceur
63      REAL flux_sparam_ind(klon), flux_sparam_bb(klon)
64c========================= LOCAL VARIABLES =============================
65      INTEGER i, k, kkk_cont(klon), kkk_expl(klon)
66      REAL zalt(klon,klev), zaltmid(klon,klev)
67      REAL zzdz
68c------------------------- SULFUR emissions ----------------------------
69      REAL lmt_so2volc_cont(klon)  ! emissions so2 volcan (continuous)
70      REAL lmt_altvolc_cont(klon)  ! altitude  so2 volcan (continuous)
71      REAL lmt_so2volc_expl(klon)  ! emissions so2 volcan (explosive)
72      REAL lmt_altvolc_expl(klon)  ! altitude  so2 volcan (explosive)
73      REAL lmt_so2ff_l(klon)       ! emissions so2 fossil fuel (low)
74      REAL lmt_so2ff_h(klon)       ! emissions so2 fossil fuel (high)
75      REAL lmt_so2nff(klon)        ! emissions so2 non-fossil fuel
76      REAL lmt_so2bb_l(klon)       ! emissions de so2 biomass burning (low)
77      REAL lmt_so2bb_h(klon)       ! emissions de so2 biomass burning (high)
78      REAL lmt_so2ba(klon)         ! emissions de so2 bateau
79      REAL lmt_dms(klon)           ! emissions de dms
80      REAL lmt_dmsconc(klon)       ! concentration de dms oceanique
81      REAL lmt_dmsbio(klon)        ! emissions de dms bio
82      REAL lmt_h2sbio(klon)        ! emissions de h2s bio
83                       
84      EXTERNAL condsurfs, liss, nightingale
85c=========================================================================
86c Modifications introduced by NHL
87c -Variables to save fluxes were introduced
88c -lmt_so2ba was multiplied by fracso2emis in line 117
89c -scale_param_bb was introduced in line 105
90c The last two modifications were errors existing in the original version
91c=========================================================================
92c=========================================================================
93c                        LOW LEVEL EMISSIONS
94c=========================================================================
95                       
96         CALL nightingale(u_seri, v_seri, u10m_ec, v10m_ec, paprs,
97     .                 pplay, cdragh, cdragm, t_seri, q_seri, ftsol,
98     .                 tsol, pctsrf, lmt_dmsconc, lmt_dms)
99
100      IF (.not.bateau) THEN
101        DO i=1, klon     
102          lmt_so2ba(i)=0.0
103        ENDDO
104      ENDIF
105
106      DO i=1, klon
107         IF (iregion_ind(i).GT.0) THEN
108       IF(id_prec>0) source_tr(i,id_prec)=source_tr(i,id_prec)
109     .             + fracso2emis
110     .             *scale_param_ind(iregion_ind(i))*lmt_so2ff_l(i)*1.e4
111     .             +scale_param_ind(iregion_ind(i))*lmt_so2ff_l(i)*1.e4
112     .                   *frach2sofso2            ! molec/m2/s
113c
114      IF(id_fine>0) source_tr(i,id_fine)=
115     .                source_tr(i,id_fine)+(1-fracso2emis)
116     .                *scale_param_ind(iregion_ind(i))*lmt_so2ff_l(i)
117     .                *1.e4*masse_ammsulfate/RNAVO  ! g/m2/s
118c
119       IF(id_prec>0)   flux_tr(i,id_prec)=flux_tr(i,id_prec) + (
120     .               scale_param_ind(iregion_ind(i))*(lmt_so2ff_l(i)+
121     .                lmt_so2ff_h(i))
122     .                *frach2sofso2
123     .               +scale_param_ind(iregion_ind(i))*(lmt_so2ff_l(i)+
124     .                lmt_so2ff_h(i))
125     .                *fracso2emis
126     .                )*1.e4/RNAVO*masse_s*1.e3          ! mgS/m2/s
127c
128      IF(id_fine>0)  flux_tr(i,id_fine)=
129     . flux_tr(i,id_fine)+(1-fracso2emis)
130     .               *scale_param_ind(iregion_ind(i))*(lmt_so2ff_l(i)+
131     .                lmt_so2ff_h(i))
132     .               *1.e4/RNAVO*masse_ammsulfate*1.e3    ! mgS/m2/s
133c
134      flux_sparam_ind(i)=flux_sparam_ind(i)+ (1-fracso2emis)
135     .               *scale_param_ind(iregion_ind(i))*(lmt_so2ff_l(i)+
136     .                lmt_so2ff_h(i))
137     .               *1.e4/RNAVO*masse_ammsulfate*1.e3    ! mgS/m2/s
138         ENDIF
139         IF (iregion_bb(i).GT.0) THEN
140      IF(id_prec>0) source_tr(i,id_prec)=
141     .                  source_tr(i,id_prec) + fracso2emis
142     .                 *scale_param_bb(iregion_bb(i))*lmt_so2bb_l(i)
143     .                   *(1.-pctsrf(i,is_oce))*1.e4
144c
145      IF(id_fine>0)     source_tr(i,id_fine)=
146     .                   source_tr(i,id_fine)+(1-fracso2emis)
147     .                  *scale_param_bb(iregion_bb(i))*lmt_so2bb_l(i)*
148     .                   (1.-pctsrf(i,is_oce))*1.e4*
149     .                   masse_ammsulfate/RNAVO  ! g/m2/s
150c
151      IF(id_prec>0)     flux_tr(i,id_prec)=flux_tr(i,id_prec) +
152     .               (scale_param_bb(iregion_bb(i))*lmt_so2bb_l(i)
153     .                 +scale_param_bb(iregion_bb(i))*lmt_so2bb_h(i))
154     .                      * (1.-pctsrf(i,is_oce))*fracso2emis
155     .                 *1.e4/RNAVO*masse_s*1.e3          ! mgS/m2/s
156c
157      IF(id_fine>0) flux_tr(i,id_fine)=
158     .                flux_tr(i,id_fine)+(1-fracso2emis)
159     .               *(scale_param_bb(iregion_bb(i))*lmt_so2bb_l(i)
160     .                +scale_param_bb(iregion_bb(i))*lmt_so2bb_h(i))
161     .                *(1.-pctsrf(i,is_oce))
162     .                *1.e4/RNAVO*masse_ammsulfate*1.e3    ! mgS/m2/s
163c
164           flux_sparam_bb(i)=
165     .                scale_param_bb(iregion_bb(i))*(lmt_so2bb_l(i)+
166     .                                        lmt_so2bb_h(i))
167     .                      * (1.-pctsrf(i,is_oce))*fracso2emis
168     .              *1.e4/RNAVO*masse_s*1.e3          ! mgS/m2/s
169           flux_sparam_bb(i)= flux_sparam_bb(i) + (1-fracso2emis) *
170     .               (scale_param_bb(iregion_bb(i))*lmt_so2bb_l(i)+
171     .                scale_param_bb(iregion_bb(i))*lmt_so2bb_h(i))
172     .                *(1.-pctsrf(i,is_oce))
173     .               *1.e4/RNAVO*masse_ammsulfate*1.e3    ! mgS/m2/s
174         ENDIF
175      IF(id_prec>0)   source_tr(i,id_prec)=source_tr(i,id_prec)
176     .                 + fracso2emis
177     .                 *(lmt_so2ba(i)+lmt_so2nff(i))*1.e4
178     .                 +(lmt_h2sbio(i)
179     .                   +lmt_dms(i)+lmt_dmsbio(i))*1.e4            ! molec/m2/s
180c
181      IF(id_fine>0)   source_tr(i,id_fine)=source_tr(i,id_fine)
182     .                +(1-fracso2emis)
183     .                *(lmt_so2ba(i)+lmt_so2nff(i))*1.e4*
184     .                   masse_ammsulfate/RNAVO  ! g/m2/s
185c
186      IF(id_prec>0)   flux_tr(i,id_prec)=flux_tr(i,id_prec)
187     .               + (lmt_h2sbio(i)
188     .               +lmt_so2volc_cont(i)+lmt_so2volc_expl(i)
189     .                 +(lmt_so2ba(i)+lmt_so2nff(i))*fracso2emis
190     .                 +lmt_dms(i)+lmt_dmsbio(i) )
191     .              *1.e4/RNAVO*masse_s*1.e3          ! mgS/m2/s
192c
193      IF(id_fine>0)   flux_tr(i,id_fine)=flux_tr(i,id_fine)
194     .               +(1-fracso2emis)
195     .               *(lmt_so2ba(i) + lmt_so2nff(i))
196     .               *1.e4/RNAVO*masse_ammsulfate*1.e3    ! mgS/m2/s
197c
198         flux_sparam_ind(i)=flux_sparam_ind(i)+ (1-fracso2emis)
199     .               *lmt_so2nff(i)
200     .               *1.e4/RNAVO*masse_ammsulfate*1.e3    ! mgS/m2/s
201c
202      ENDDO
203
204c========================================================================
205c                        HIGH LEVEL EMISSIONS
206c========================================================================
207c  Source de SO2 volcaniques
208      DO i = 1, klon
209        kkk_cont(i)=1
210        kkk_expl(i)=1
211      ENDDO
212      DO k=1, klev-1
213      DO i = 1, klon
214        zaltmid(i,k)=zalt(i,k)+zdz(i,k)/2.
215        IF (zalt(i,k+1).LT.lmt_altvolc_cont(i)) kkk_cont(i)=k+1
216        IF (zalt(i,k+1).LT.lmt_altvolc_expl(i)) kkk_expl(i)=k+1
217      ENDDO
218      ENDDO
219      IF(id_prec>0) THEN
220      DO i = 1, klon
221        tr_seri(i,kkk_cont(i),id_prec)=tr_seri(i,kkk_cont(i),id_prec) +
222     .               lmt_so2volc_cont(i)/zdz(i,kkk_cont(i))/100.*pdtphys
223        tr_seri(i,kkk_expl(i),id_prec)=tr_seri(i,kkk_expl(i),id_prec) +
224     .               lmt_so2volc_expl(i)/zdz(i,kkk_expl(i))/100.*pdtphys
225      ENDDO
226      ENDIF                                           
227c  Sources hautes de SO2     
228     
229c
230c--only GEIA SO2 emissions has high emissions
231c--unit: molec/cm2/s divided by layer height (in cm) multiplied by timestep
232c
233      k=2                             !introducing emissions in level 2
234      DO i = 1, klon
235c
236         IF (iregion_bb(i).GT.0) THEN
237      IF(id_prec>0)   tr_seri(i,k,id_prec)=
238     .              tr_seri(i,k,id_prec) + fracso2emis
239     .              *scale_param_bb(iregion_bb(i))*lmt_so2bb_h(i)
240     .              /zdz(i,k)/100.*pdtphys
241c
242      IF(id_fine>0)     tr_seri(i,k,id_fine)=tr_seri(i,k,id_fine)
243     .              + (1.-fracso2emis)
244     .              *scale_param_bb(iregion_bb(i))*lmt_so2bb_h(i)
245     .              *masse_ammsulfate/RNAVO/zdz(i,k)/100.*pdtphys   !g/cm3
246         ENDIF
247         IF (iregion_ind(i).GT.0) THEN
248       IF(id_prec>0)  tr_seri(i,k,id_prec)=
249     .              tr_seri(i,k,id_prec) + (fracso2emis
250     .              *scale_param_ind(iregion_ind(i))*lmt_so2ff_h(i)
251     .              + frach2sofso2
252     .              *scale_param_ind(iregion_ind(i))*lmt_so2ff_h(i))
253     .              /zdz(i,k)/100.*pdtphys
254c
255       IF(id_fine>0)    tr_seri(i,k,id_fine)=tr_seri(i,k,id_fine)
256     .               + (1.-fracso2emis)
257     .              *scale_param_ind(iregion_ind(i))*lmt_so2ff_h(i)
258     .              *masse_ammsulfate/RNAVO/zdz(i,k)/100.*pdtphys   !g/cm3
259         ENDIF
260c
261      ENDDO
262
263      END
Note: See TracBrowser for help on using the repository browser.