1 | subroutine hazesource(ngrid,nlayer,nq,ptimestep, & |
---|
2 | pplev,flusurf,mu0,zdq_source) |
---|
3 | |
---|
4 | use radinc_h, only : naerkind |
---|
5 | use comgeomfi_h |
---|
6 | implicit none |
---|
7 | |
---|
8 | !================================================================== |
---|
9 | ! Purpose |
---|
10 | ! ------- |
---|
11 | ! Surface source of haze particle |
---|
12 | ! |
---|
13 | ! Inputs |
---|
14 | ! ------ |
---|
15 | ! ngrid Number of vertical columns |
---|
16 | ! nlayer Number of layers |
---|
17 | ! pplay(ngrid,nlayer) Pressure layers |
---|
18 | ! pplev(ngrid,nlayer+1) Pressure levels |
---|
19 | ! |
---|
20 | ! Outputs |
---|
21 | ! ------- |
---|
22 | ! zdq_source |
---|
23 | ! |
---|
24 | ! Authors |
---|
25 | ! ------- |
---|
26 | ! Tanguy Bertrand |
---|
27 | ! |
---|
28 | !================================================================== |
---|
29 | |
---|
30 | #include "dimensions.h" |
---|
31 | #include "dimphys.h" |
---|
32 | #include "comcstfi.h" |
---|
33 | #include "surfdat.h" |
---|
34 | #include "comvert.h" |
---|
35 | #include "callkeys.h" |
---|
36 | #include "tracer.h" |
---|
37 | |
---|
38 | !----------------------------------------------------------------------- |
---|
39 | ! Arguments |
---|
40 | |
---|
41 | INTEGER,INTENT(IN) :: ngrid, nlayer, nq |
---|
42 | REAL,INTENT(IN) :: ptimestep |
---|
43 | REAL,INTENT(IN) :: pplev(ngrid,nlayer+1) |
---|
44 | REAL,INTENT(IN) :: flusurf(ngrid,nq) ! flux cond/sub kg.m-2.s-1 |
---|
45 | REAL,INTENT(IN) :: mu0(ngrid) ! cosinus of solar incident flux |
---|
46 | REAL,INTENT(OUT) :: zdq_source(ngrid,nlayer,nq) |
---|
47 | !----------------------------------------------------------------------- |
---|
48 | ! Local variables |
---|
49 | |
---|
50 | INTEGER l,ig,iq |
---|
51 | CHARACTER(len=10) :: tname |
---|
52 | real massfix(ngrid) ! |
---|
53 | real tholflux(ngrid) ! |
---|
54 | !----------------------------------------------------------------------- |
---|
55 | |
---|
56 | !! From callphys.def : kfix, fracsource, latsource, mode_hs |
---|
57 | |
---|
58 | !! 1) Get mass of the layer (kg/m2) where haze is injected |
---|
59 | massfix(:)=0.0 |
---|
60 | zdq_source(:,:,:)=0.0 |
---|
61 | DO ig=1,ngrid |
---|
62 | DO l= 1,kfix |
---|
63 | massfix(ig)= massfix(ig) + (pplev(ig,l)-pplev(ig,l+1))/g |
---|
64 | ENDDO |
---|
65 | ENDDO |
---|
66 | |
---|
67 | |
---|
68 | SELECT CASE (mode_hs) ! mode haze source |
---|
69 | !!! Source haze in SP: 0.02 pourcent when n2 sublimes |
---|
70 | CASE (0) |
---|
71 | DO ig=1,ngrid |
---|
72 | IF (phisfi(ig).le.-1500.) THEN ! SP |
---|
73 | IF (lati(ig)*180./pi.ge.latsource.and.lati(ig)*180./pi.lt.65.) THEN |
---|
74 | |
---|
75 | IF (flusurf(ig,igcm_n2).lt.0) THEN |
---|
76 | IF (mu0(ig).gt.0.6) THEN |
---|
77 | !! flux tholins kg/m2/s |
---|
78 | tholflux(ig)=-flusurf(ig,igcm_n2)*fracsource |
---|
79 | DO iq=1,nq |
---|
80 | tname=noms(iq) |
---|
81 | if (tname(1:4).eq."haze") then |
---|
82 | DO l= 1,kfix |
---|
83 | zdq_source(ig,l,iq)=tholflux(ig)/massfix(ig) |
---|
84 | ENDDO |
---|
85 | endif |
---|
86 | ENDDO |
---|
87 | ENDIF ! mu0 |
---|
88 | ENDIF ! flusurf |
---|
89 | |
---|
90 | ENDIF |
---|
91 | ENDIF |
---|
92 | ENDDO |
---|
93 | |
---|
94 | !!! Source haze at the BTD |
---|
95 | CASE (1) |
---|
96 | DO ig=1,ngrid |
---|
97 | IF (phisfi(ig).ge.0.) THEN |
---|
98 | !print*,'lat/lon=',lati(ig)*180./pi,long(ig)*180./pi |
---|
99 | IF (lati(ig)*180./pi.ge.0..and.lati(ig)*180./pi.lt.5.) THEN |
---|
100 | IF (long(ig)*180./pi.ge.-50..and.long(ig)*180./pi.lt.-45.) THEN |
---|
101 | !print*,'YES',flusurf(ig,igcm_ch4_ice) |
---|
102 | |
---|
103 | IF (flusurf(ig,igcm_ch4_ice).lt.0) THEN |
---|
104 | !! flux tholins kg/m2/s |
---|
105 | tholflux(ig)=-flusurf(ig,igcm_ch4_ice)*fracsource |
---|
106 | DO iq=1,nq |
---|
107 | tname=noms(iq) |
---|
108 | if (tname(1:4).eq."haze") then |
---|
109 | DO l= 1,kfix |
---|
110 | zdq_source(ig,l,iq)=tholflux(ig)/massfix(ig) |
---|
111 | ENDDO |
---|
112 | endif |
---|
113 | ENDDO |
---|
114 | ENDIF ! flusurf |
---|
115 | |
---|
116 | ENDIF |
---|
117 | ENDIF |
---|
118 | ENDIF |
---|
119 | ENDDO |
---|
120 | |
---|
121 | !!! Basic source(geysers) |
---|
122 | CASE (2) |
---|
123 | DO ig=1,ngrid |
---|
124 | IF ((abs(lati(ig)*180./pi-latsource).le.abs(lati(1)-lati(2))*180./pi) & |
---|
125 | .and.(abs(long(ig)*180./pi-lonsource).le.abs(long(3)-long(4))*180./pi)) THEN |
---|
126 | |
---|
127 | !! flux tholins kg/m2/s |
---|
128 | tholflux(ig)=fracsource |
---|
129 | DO iq=1,nq |
---|
130 | tname=noms(iq) |
---|
131 | if (tname(1:4).eq."haze") then |
---|
132 | DO l= 1,kfix |
---|
133 | zdq_source(ig,l,iq)=tholflux(ig)/massfix(ig) |
---|
134 | ENDDO |
---|
135 | endif |
---|
136 | ENDDO |
---|
137 | ENDIF |
---|
138 | ENDDO |
---|
139 | |
---|
140 | |
---|
141 | CASE DEFAULT |
---|
142 | write(*,*) 'STOP in hazesource.F90: mod_hs not found' |
---|
143 | stop |
---|
144 | END SELECT |
---|
145 | |
---|
146 | |
---|
147 | end subroutine hazesource |
---|
148 | |
---|