source: trunk/LMDZ.GENERIC/libf/phystd/interpolateH2He.F90 @ 3562

Last change on this file since 3562 was 3310, checked in by sglmd, 8 months ago

include updated H2-H2 and H2-He CIA files by G. Milcareck for equilibrium case

File size: 8.3 KB
RevLine 
[878]1     subroutine interpolateH2He(wn,temp,presH2,presHe,abcoef,firstcall,ind)
[716]2
3!==================================================================
4!     
5!     Purpose
6!     -------
7!     Calculates the H2-He CIA opacity, using a lookup table from
8!     HITRAN (2011)
9!
10!     Authors
11!     -------
12!     R. Wordsworth (2011)
13!     
14!==================================================================
15
[2655]16      use callkeys_mod, only: H2orthopara_mixture, strictboundcia
[716]17      use datafile_mod, only: datadir
[2662]18      use mod_phys_lmdz_para, only : is_master
[873]19
[716]20      implicit none
21
22      ! input
23      double precision wn                 ! wavenumber             (cm^-1)
24      double precision temp               ! temperature            (Kelvin)
25      double precision presH2             ! H2 partial pressure    (Pascals)
26      double precision presHe             ! He partial pressure    (Pascals)
27
28      ! output
29      double precision abcoef             ! absorption coefficient (m^-1)
30
31      integer nS,nT
32      parameter(nT=10)
33
34      double precision, parameter :: losch = 2.6867774e19
35      ! Loschmit's number (molecule cm^-3 at STP)
36      ! converts cm^5 molecule^-2 --> cm^-1 amagat^-2
37      ! see Richard et al. 2011, JQSRT for details
38
39      double precision amagatH2
40      double precision amagatHe
41      double precision temp_arr(nT)
[2667]42      double precision, dimension(:),   allocatable :: wn_arr
43      double precision, dimension(:,:), allocatable :: abs_arr
[716]44
45      integer k,iT
46      logical firstcall
47
[3310]48      save ns, wn_arr, temp_arr, abs_arr !read by master
[716]49
50      character*100 dt_file
51      integer strlen,ios
52
53      character(LEN=*), parameter :: fmat1 = "(A20,F10.3,F10.3,I7,F7.1,E10.3,F5.3)"
54
55      character*20 bleh
[2870]56      double precision blah, Ttemp, ztemp
57
[716]58      integer nres
59
[878]60      integer ind
[716]61
[2870]62      ztemp=temp
[2667]63
64      if (H2orthopara_mixture .eq. "hot") then ! .and. (temp .gt. 3000.0)) then
65        ! print*,"We're in the Hot Jupiter case "
[2870]66        if (ztemp .gt. 9900) then
[2667]67          if (strictboundcia) then
68            if (is_master) then
69              print*,'Your temperatures are too high for this H2-He CIA dataset (Hot Jupiter case).'
70              print*,'Please run mixed H2-He atmospheres below T = 9900.0 K.'     
71            endif !is_master
72            stop
73          else
[2870]74            !if (is_master) then
75            !  print*,'Your temperatures are too high for this H2-He CIA dataset (Hot Jupiter case)'
76            !  print*,'you have chosen strictboundcia = ', strictboundcia
77            !  print*,'*********************************************************'
78            !  print*,' we allow model to continue but with temp = 9900.0 K     '
79            !  print*,'  ...       for H2-He CIA dataset          ...           '
80            !  print*,'  ... we assume we know what you are doing ...           '
81            !  print*,'*********************************************************'
82            !endif !is_master
83            ztemp = 9900.0
[2667]84          endif ! of stricbound cia
85        endif ! temp .gt. 9900
86      else !not in Hot Jupiter case
[2870]87        if(ztemp.gt.400)then
[2667]88          if (strictboundcia) then
89            if (is_master) then
90              print*,'Your temperatures are too high for this H2-He CIA dataset.'
91              print*,'Please run mixed H2-He atmospheres below T = 400 K.'     
92            endif ! is_master
93            stop
94          else
[2870]95            !if (is_master) then
96            !  print*,'Your temperatures are too high for this H2-He CIA dataset'
97            !  print*,'you have chosen strictboundcia = ', strictboundcia
98            !  print*,'*********************************************************'
99            !  print*,' we allow model to continue but with temp = 400          '
100            !  print*,'  ...       for H2-He CIA dataset          ...           '
101            !  print*,'  ... we assume we know what you are doing ...           '
102            !  print*,'*********************************************************'
103            !endif ! is_master
104            ztemp = 400
[2667]105          endif !of strictboundcia
[2870]106        elseif(ztemp.lt.40)then
[2667]107          if (strictboundcia) then
108            if (is_master) then
109              print*,'Your temperatures are too low for this H2-He CIA dataset.'
110              print*,'Please run mixed H2-He atmospheres above T = 40 K.' 
111            endif ! is_master   
112            stop
113          else
[2870]114            !if (is_master) then
115            !  print*,'Your temperatures are too low for this H2-He CIA dataset'
116            !  print*,'you have chosen strictboundcia = ', strictboundcia
117            !  print*,'*********************************************************'
118            !  print*,' we allow model to continue but with temp = 40           '
119            !  print*,'  ...       for H2-He CIA dataset          ...           '
120            !  print*,'  ... we assume we know what you are doing ...           '
121            !  print*,'*********************************************************'
122            !endif ! is_master
123            ztemp = 40
[2667]124          endif !of strictboundcia
125        endif ! of (temp .gt. 400)
126      endif ! of (H2orthopara_mixture .eq. "hot")
127
128
[716]129      amagatH2 = (273.15/temp)*(presH2/101325.0)
130      amagatHe = (273.15/temp)*(presHe/101325.0)
131
132      if(firstcall)then ! called by sugas_corrk only
[2662]133         if (is_master) print*,'----------------------------------------------------'
134         if (is_master) print*,'Initialising H2-He continuum from HITRAN database...'
[716]135
136!     1.1 Open the ASCII files
[2655]137         if (H2orthopara_mixture.eq."normal") then
[3310]138            nS = 2428
139            dt_file=TRIM(datadir)//'/continuum_data/H2-He_norm_2011.cia'
[2655]140         else if (H2orthopara_mixture.eq."equilibrium") then
[3310]141            nS = 1501
142            dt_file=TRIM(datadir)//'/continuum_data/H2-He_eq_0-15000cm-1_40-400K.cia'
[2667]143         else if (H2orthopara_mixture .eq. "hot") then !we use a dataset than can go as high as 9900 K
[3310]144            nS = 19981
145            dt_file=TRIM(datadir)//'/continuum_data/H2-He_2011.cia'
[2655]146         endif
[3310]147
148         if (.not. allocated(wn_arr)) allocate(wn_arr(nS))
149         if (.not. allocated(abs_arr)) allocate(abs_arr(nS,nT))
150
[1315]151         
152!$OMP MASTER
[716]153         open(33,file=dt_file,form='formatted',status='old',iostat=ios)
154         if (ios.ne.0) then        ! file not found
[2662]155           if (is_master) then
156             write(*,*) 'Error from interpolateH2He'
157             write(*,*) 'Data file ',trim(dt_file),' not found.'
158             write(*,*) 'Check that your path to datagcm:',trim(datadir)
159             write(*,*) 'is correct. You can change it in callphys.def with:'
160             write(*,*) 'datadir = /absolute/path/to/datagcm'
161             write(*,*) 'Also check that the continuum data is there.'
162           endif
[716]163           call abort
164         else
165
166            do iT=1,nT
167
168               read(33,fmat1) bleh,blah,blah,nres,Ttemp
169               if(nS.ne.nres)then
[2662]170                  if (is_master) then
171                    print*,'Resolution given in file: ',trim(dt_file)
172                    print*,'is ',nres,', which does not match nS.'
173                    print*,'Please adjust nS value in interpolateH2He.F90'
174                  endif
[716]175                  stop
176               endif
177               temp_arr(iT)=Ttemp
178
179               do k=1,nS
180                  read(33,*) wn_arr(k),abs_arr(k,it)
181               end do
182
183            end do
184     
185         endif
186         close(33)
[1315]187!$OMP END MASTER
188!$OMP BARRIER
[2662]189         if (is_master) then
190           print*,'interpolateH2He: At wavenumber ',wn,' cm^-1'
191           print*,'   temperature                 ',temp,' K'
192           print*,'   H2 partial pressure         ',presH2,' Pa'
193           print*,'   and He partial pressure     ',presHe,' Pa'
194         endif
[3310]195      endif !! if (firstcall)
[716]196
[2870]197         call bilinearbig(nS,nT,wn_arr,temp_arr,abs_arr,wn,ztemp,abcoef,ind)
[716]198
[873]199         !print*,'the absorption is ',abcoef,' cm^5 molecule^-2'
200         !print*,'or ',abcoef*losch**2,' cm^-1 amagat^-2'
201
[716]202         abcoef=abcoef*losch**2*100.0*amagatH2*amagatHe ! convert to m^-1
203
[873]204         !print*,'We have ',amagatH2,' amagats of H2'
205         !print*,'and     ',amagatHe,' amagats of He'
206         !print*,'So the absorption is ',abcoef,' m^-1'
[716]207
208         ! unlike for Rayleigh scattering, we do not currently weight by the BB function
209         ! however our bands are normally thin, so this is no big deal.
210
211      return
212    end subroutine interpolateH2He
Note: See TracBrowser for help on using the repository browser.