source: trunk/LMDZ.VENUS/libf/phyvenus/compo_hedin_mod2.F90 @ 3094

Last change on this file since 3094 was 2686, checked in by slebonnois, 3 years ago

SL+AM : various corrections prior to VCD 2.0

File size: 5.2 KB
Line 
1  module compo_hedin83_mod2
2  !!05/2013 Laura Salmi
3  !!03/2014 revision Gabriella Gilli
4  !!Calcul des vmr pour CO2, CO, O, N et N2 en s'appuyant sur les tables donnees dans l'article Hedin (1983)
5  implicit none     
6  contains
7
8        subroutine compo_hedin83_init2
9
10        implicit none
11#include "YOMCST.h"
12#include "hedin.h"
13
14        REAL :: alpha
15        REAL :: mu(musize),z(zsize)
16        REAL :: T(musize,zsize)
17        REAL :: ntot(musize, zsize)
18        REAL :: nco2(musize, zsize),nn2(musize, zsize)
19        REAL :: nco(musize, zsize), no(musize, zsize)
20        REAL :: nn(musize, zsize),nhe(musize, zsize)
21        REAL :: sig_co2(zsize), sig_o(zsize), sig_n2(zsize),sig_co(zsize)
22        integer :: i,j
23        REAL factor_ox3p
24
25
26! Initialisation des sza et de l'atlitude
27        mu(1)=1.
28        do i=2,musize-1
29              mu(i)=mu(i-1)-(2./(1.*musize-1.))
30        enddo
31        mu(musize)=-1.
32
33        z(1)=100.
34        do i=2,zsize
35              z(i)=z(i-1)+5.
36        enddo
37
38!Lecture des variables (tables de Hedin (1983))
39
40       open( 11, file='HIGHATM/noon.txt')
41      read (11, *)
42      do i=1,zsize
43         read (11,*) T(1,i), nco2(1,i),no(1,i), nco(1,i), nn2(1,i),nn(1,i), nhe(1,i)
44      enddo
45      close (11)
46       open( 12, file='HIGHATM/midnight.txt')
47      read (12, *)
48      do i=1,zsize
49         read (12,*) T(musize,i), nco2(musize,i),no(musize,i), nco(musize,i), nn2(musize,i), nn(musize,i), nhe(musize,i)
50      enddo
51      close (12)
52
53!Dependance en SZA
54
55        do i=2,musize-1
56          alpha=0.5*tanh(4.*(mu(i)-0.5))+0.5
57          do j=1,zsize
58                T(i,j)=alpha*T(1,j)+(1.-alpha)*T(musize,j)
59                nco2(i,j)=alpha*nco2(1,j)+(1.-alpha)*nco2(musize,j)
60                no(i,j)=alpha*no(1,j)+(1.-alpha)*no(musize,j)
61                nco(i,j)=alpha*nco(1,j)+(1.-alpha)*nco(musize,j)
62                nn2(i,j)=alpha*nn2(1,j)+(1.-alpha)*nn2(musize,j)
63                nn(i,j)=alpha*nn(1,j)+(1.-alpha)*nn(musize,j)
64                nhe(i,j)=alpha*nhe(1,j)+(1.-alpha)*nhe(musize,j)
65
66          enddo
67         enddo
68
69!! Test: effect of varying atomic oxygen abundances
70
71         factor_ox3p = 2.
72           DO i=1,musize
73             DO j=1,zsize
74               no(i,j)    =  factor_ox3p* no(i,j)
75             ENDDO
76           ENDDO
77
78
79
80
81
82
83!! Conversion en volume mixture ratio
84
85           do i=1,musize
86                  do j=1,zsize
87                  ntot(i,j)=nco2(i,j)+nn2(i,j)+nn(i,j)+nco(i,j)+no(i,j)+nhe(i,j)+2e-3*nco2(i,j) !for o2 and no si on change la proportion changer aussi dans euvheat
88                  pres_hedin(i,j)=ntot(i,j)*1e6*RKBOL*T(i,j)
89                 co2_hedin(i,j)=nco2(i,j)/ntot(i,j)
90                 co_hedin(i,j)=nco(i,j)/ntot(i,j)
91                 n2_hedin(i,j)=nn2(i,j)/ntot(i,j)
92                 n_hedin(i,j)=nn(i,j)/ntot(i,j)
93                 o_hedin(i,j)=no(i,j)/ntot(i,j)
94                 mu_hedin(i)=mu(i)
95
96          enddo
97         enddo
98
99!           print*, pres_hedin(1,:)
100!           print*, ' '
101!           print*, T(musize,:) 
102!           print*, ' '
103!            print*, T(1,:)
104!               print*, nn2(10,:)
105!               print*, ' '
106!               print*,nco2(10,:)
107!               print*, ' '
108!               print*,nco(10,:)
109!               print*, ' '
110!               print*,no(10,:)
111
112!               print*, " "
113!               print*, ntot(1,:)
114!               print*, co2_hedin(1,:)
115!               print*, mu_hedin(:)
116!               print*, " "
117!               print*, n2_hedin(1,:)
118!               print*, " "
119!               print*, o_hedin(1,:)
120!               print*, " "
121!               print*, co_hedin(1,:)
122!               print*, ' '
123!           stop
124
125       end subroutine compo_hedin83_init2
126
127
128
129SUBROUTINE compo_hedin83_mod(pression,mu0,co2vmr_gcm,covmr_gcm,ovmr_gcm,n2vmr_gcm,nvmr_gcm)
130
131
132!!Interpolation des profils de Hedin (1983) sur la grille du modele
133
134        use dimphy
135        implicit none
136        include 'hedin.h'
137        REAL, intent(in) :: pression(klon,klev), mu0(klon)
138        integer :: i,k,iz,jmu,jz,z_ok, mu_ok
139        REAL :: factp(klev)
140        REAL, intent(out) :: co2vmr_gcm(klon,klev),covmr_gcm(klon,klev)
141        REAL, intent(out) ::ovmr_gcm(klon,klev), n2vmr_gcm(klon,klev), nvmr_gcm(klon,klev)
142        REAL :: ang0, ang_hedin
143!        print*, pres_hedin(1,:)
144!       print*, ' '
145!       print*, pres_hedin(10,:)
146!       stop
147       do i=1,klon
148         ang0=acos(mu0(i))*180./3.1415
149         do jmu=1,musize
150            ang_hedin=acos(mu_hedin(jmu))*180./3.1415
151         if (ang_hedin.le.ang0 ) then   
152!          if (mu_hedin(jmu) .le. mu0(i) ) then 
153              mu_ok = jmu         
154         endif
155         enddo
156!       print*, ang_hedin, ang0, mu_ok
157!       STOP
158
159        do k=1,klev
160         z_ok=2         
161         do jz=2,zsize-1
162           if (pres_hedin(mu_ok,jz).ge.pression(i,k)) then
163                 z_ok = jz+1
164           endif     
165          enddo
166
167        if (pres_hedin(mu_ok,zsize).ge.pression(i,k)) then
168! avoid the extrapolation... Problem when p < 1.E-6 Pa...
169         factp(k) = 1.
170        else
171         factp(k) = (log10(pression(i,k))-log10(pres_hedin(mu_ok,z_ok-1)))/(log10(pres_hedin(mu_ok,z_ok))-log10(pres_hedin(mu_ok,z_ok-1)))
172        endif
173
174        ovmr_gcm(i,k) = (10.**(log10(o_hedin(mu_ok,z_ok))*factp(k)+log10(o_hedin(mu_ok,z_ok-1))*(1.-factp(k))))
175        n2vmr_gcm(i,k) = 10.**(log10(n2_hedin(mu_ok,z_ok))*factp(k)+log10(n2_hedin(mu_ok,z_ok-1))*(1.-factp(k)))
176        nvmr_gcm(i,k) = 10.**(log10(n_hedin(mu_ok,z_ok))*factp(k)+log10(n_hedin(mu_ok,z_ok-1))*(1.-factp(k)))
177        covmr_gcm(i,k) = 10.**(log10(co_hedin(mu_ok,z_ok))*factp(k)+log10(co_hedin(mu_ok,z_ok-1))*(1.-factp(k)))
178        co2vmr_gcm(i,k) =1.- (ovmr_gcm(i,k)+n2vmr_gcm(i,k)+covmr_gcm(i,k))
179!        co2vmr_gcm(i,k) =1.- (ovmr_gcm(i,k)+n2vmr_gcm(i,k)+covmr_gcm(i,k)+nvmr_gcm(i,k))
180
181        enddo   
182       enddo
183!      stop
184       
185
186END SUBROUTINE compo_hedin83_mod
187
188
189  end module compo_hedin83_mod2
190
Note: See TracBrowser for help on using the repository browser.