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

Last change on this file since 2937 was 2870, checked in by jleconte, 2 years ago

Removed too many prints in CIA interpolation + correction for hot H2H2

File size: 8.2 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
[1315]48      save 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 "
66        nS = 19981
[2870]67        if (ztemp .gt. 9900) then
[2667]68          if (strictboundcia) then
69            if (is_master) then
70              print*,'Your temperatures are too high for this H2-He CIA dataset (Hot Jupiter case).'
71              print*,'Please run mixed H2-He atmospheres below T = 9900.0 K.'     
72            endif !is_master
73            stop
74          else
[2870]75            !if (is_master) then
76            !  print*,'Your temperatures are too high for this H2-He CIA dataset (Hot Jupiter case)'
77            !  print*,'you have chosen strictboundcia = ', strictboundcia
78            !  print*,'*********************************************************'
79            !  print*,' we allow model to continue but with temp = 9900.0 K     '
80            !  print*,'  ...       for H2-He CIA dataset          ...           '
81            !  print*,'  ... we assume we know what you are doing ...           '
82            !  print*,'*********************************************************'
83            !endif !is_master
84            ztemp = 9900.0
[2667]85          endif ! of stricbound cia
86        endif ! temp .gt. 9900
87      else !not in Hot Jupiter case
88        nS = 2428
[2870]89        if(ztemp.gt.400)then
[2667]90          if (strictboundcia) then
91            if (is_master) then
92              print*,'Your temperatures are too high for this H2-He CIA dataset.'
93              print*,'Please run mixed H2-He atmospheres below T = 400 K.'     
94            endif ! is_master
95            stop
96          else
[2870]97            !if (is_master) then
98            !  print*,'Your temperatures are too high for this H2-He CIA dataset'
99            !  print*,'you have chosen strictboundcia = ', strictboundcia
100            !  print*,'*********************************************************'
101            !  print*,' we allow model to continue but with temp = 400          '
102            !  print*,'  ...       for H2-He CIA dataset          ...           '
103            !  print*,'  ... we assume we know what you are doing ...           '
104            !  print*,'*********************************************************'
105            !endif ! is_master
106            ztemp = 400
[2667]107          endif !of strictboundcia
[2870]108        elseif(ztemp.lt.40)then
[2667]109          if (strictboundcia) then
110            if (is_master) then
111              print*,'Your temperatures are too low for this H2-He CIA dataset.'
112              print*,'Please run mixed H2-He atmospheres above T = 40 K.' 
113            endif ! is_master   
114            stop
115          else
[2870]116            !if (is_master) then
117            !  print*,'Your temperatures are too low for this H2-He CIA dataset'
118            !  print*,'you have chosen strictboundcia = ', strictboundcia
119            !  print*,'*********************************************************'
120            !  print*,' we allow model to continue but with temp = 40           '
121            !  print*,'  ...       for H2-He CIA dataset          ...           '
122            !  print*,'  ... we assume we know what you are doing ...           '
123            !  print*,'*********************************************************'
124            !endif ! is_master
125            ztemp = 40
[2667]126          endif !of strictboundcia
127        endif ! of (temp .gt. 400)
128      endif ! of (H2orthopara_mixture .eq. "hot")
129
130      if (.not. allocated(wn_arr)) allocate(wn_arr(nS))
131      if (.not. allocated(abs_arr)) allocate(abs_arr(nS,nT))
132
[716]133      amagatH2 = (273.15/temp)*(presH2/101325.0)
134      amagatHe = (273.15/temp)*(presHe/101325.0)
135
136      if(firstcall)then ! called by sugas_corrk only
[2662]137         if (is_master) print*,'----------------------------------------------------'
138         if (is_master) print*,'Initialising H2-He continuum from HITRAN database...'
[716]139
140!     1.1 Open the ASCII files
[2655]141         if (H2orthopara_mixture.eq."normal") then
142           dt_file=TRIM(datadir)//'/continuum_data/H2-He_norm_2011.cia'
143         else if (H2orthopara_mixture.eq."equilibrium") then
144           dt_file=TRIM(datadir)//'/continuum_data/H2-He_eq_2011.cia'
[2667]145         else if (H2orthopara_mixture .eq. "hot") then !we use a dataset than can go as high as 9900 K
146          dt_file=TRIM(datadir)//'/continuum_data/H2-He_2011.cia'
[2655]147         endif
[1315]148         
149!$OMP MASTER
[716]150         open(33,file=dt_file,form='formatted',status='old',iostat=ios)
151         if (ios.ne.0) then        ! file not found
[2662]152           if (is_master) then
153             write(*,*) 'Error from interpolateH2He'
154             write(*,*) 'Data file ',trim(dt_file),' not found.'
155             write(*,*) 'Check that your path to datagcm:',trim(datadir)
156             write(*,*) 'is correct. You can change it in callphys.def with:'
157             write(*,*) 'datadir = /absolute/path/to/datagcm'
158             write(*,*) 'Also check that the continuum data is there.'
159           endif
[716]160           call abort
161         else
162
163            do iT=1,nT
164
165               read(33,fmat1) bleh,blah,blah,nres,Ttemp
166               if(nS.ne.nres)then
[2662]167                  if (is_master) then
168                    print*,'Resolution given in file: ',trim(dt_file)
169                    print*,'is ',nres,', which does not match nS.'
170                    print*,'Please adjust nS value in interpolateH2He.F90'
171                  endif
[716]172                  stop
173               endif
174               temp_arr(iT)=Ttemp
175
176               do k=1,nS
177                  read(33,*) wn_arr(k),abs_arr(k,it)
178               end do
179
180            end do
181     
182         endif
183         close(33)
[1315]184!$OMP END MASTER
185!$OMP BARRIER
[2662]186         if (is_master) then
187           print*,'interpolateH2He: At wavenumber ',wn,' cm^-1'
188           print*,'   temperature                 ',temp,' K'
189           print*,'   H2 partial pressure         ',presH2,' Pa'
190           print*,'   and He partial pressure     ',presHe,' Pa'
191         endif
[873]192      endif
[716]193
[2870]194         call bilinearbig(nS,nT,wn_arr,temp_arr,abs_arr,wn,ztemp,abcoef,ind)
[716]195
[873]196         !print*,'the absorption is ',abcoef,' cm^5 molecule^-2'
197         !print*,'or ',abcoef*losch**2,' cm^-1 amagat^-2'
198
[716]199         abcoef=abcoef*losch**2*100.0*amagatH2*amagatHe ! convert to m^-1
200
[873]201         !print*,'We have ',amagatH2,' amagats of H2'
202         !print*,'and     ',amagatHe,' amagats of He'
203         !print*,'So the absorption is ',abcoef,' m^-1'
[716]204
205         ! unlike for Rayleigh scattering, we do not currently weight by the BB function
206         ! however our bands are normally thin, so this is no big deal.
207
208      return
209    end subroutine interpolateH2He
Note: See TracBrowser for help on using the repository browser.