source: LMDZ6/trunk/libf/phylmd/Dust/precuremission.f90 @ 5450

Last change on this file since 5450 was 5337, checked in by Laurent Fairhead, 5 weeks ago

Getting rid of dependance to dynamics

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