source: trunk/LMDZ.GENERIC/libf/phystd/interpolateH2H2.F90 @ 2655

Last change on this file since 2655 was 2655, checked in by gmilcareck, 3 years ago

Major changes to CIA interpolation:
1) Add contribution from CH4 (H2-CH4,He-CH4,CH4-CH4) ;
2) Add some tests before interpolation for H2, He and CH4 ;
3) Add the possibility to choose between a normal or equilibrium ortho:para
fraction for CIA H2;
4) Change "strictboundH2H2cia" to the generic "strictboundcia" for H2,He,CH4.
It can be added for others CIA (N2,H,CO2...) if you want.

  • Property svn:executable set to *
File size: 7.4 KB
Line 
1
2     subroutine interpolateH2H2(wn,temp,pres,abcoef,firstcall,ind)
3
4!==================================================================
5!     
6!     Purpose
7!     -------
8!     Calculates the H2-H2 CIA opacity, using a lookup table from
9!     HITRAN (2011 or later)
10!
11!     Authors
12!     -------
13!     R. Wordsworth (2011)
14!
15!     + J.Vatant d'Ollone (2019)
16!        - Enable updated HITRAN file (Karman2019,Fletcher2018)
17!           (2018 one should be default for giant planets)
18!==================================================================
19
20      use callkeys_mod, only: versH2H2cia, strictboundcia, H2orthopara_mixture
21      use datafile_mod, only: datadir
22
23      implicit none
24
25      ! input
26      double precision wn                 ! wavenumber             (cm^-1)
27      double precision temp               ! temperature            (Kelvin)
28      double precision pres               ! pressure               (Pascals)
29
30      ! output
31      double precision abcoef             ! absorption coefficient (m^-1)
32
33      integer nS,nT
34      parameter(nT=10)
35
36      double precision, parameter :: losch = 2.6867774e19
37      ! Loschmit's number (molecule cm^-3 at STP)
38      ! converts cm^5 molecule^-2 --> cm^-1 amagat^-2
39      ! see Richard et al. 2011, JQSRT for details
40
41      double precision amagat
42      double precision temp_arr(nT)
43     
44      double precision, dimension(:),   allocatable :: wn_arr
45      double precision, dimension(:,:), allocatable :: abs_arr
46
47      integer k,iT
48      logical firstcall
49
50      save nS, wn_arr, temp_arr, abs_arr !read by master
51
52      character*100 dt_file
53      integer ios
54
55      character(LEN=*), parameter :: fmat11 = "(A20,F10.3,F10.3,I7,F7.1,E10.3,F5.3)"
56      character(LEN=*), parameter :: fmat18 = "(A12,A3,A5,F10.6,F10.4,I7,F7.3,E10.3,F5.3)"
57
58      character*20 bleh
59      double precision blah, Ttemp
60      integer nres
61
62      integer ind
63     
64      if(temp.gt.400)then
65        if (strictboundcia) then
66          print*,'Your temperatures are too high for this H2-H2 CIA dataset. If you '
67          print*,'really want to run simulations with hydrogen at T > 400 K, contact'
68          print*,'Robin Wordsworth [rwordsworth@uchicago.edu].'
69          stop
70        else
71          print*,'Your temperatures are too high for this H2-H2 CIA dataset'
72          print*,'you have chosen strictboundcia = ', strictboundcia
73          print*,'*********************************************************'
74          print*,' we allow model to continue but with temp = 400          '
75          print*,'  ...       for H2-H2 CIA dataset          ...           '
76          print*,'  ... we assume we know what you are doing ...           '
77          print*,'*********************************************************'
78          temp = 400
79        endif
80      elseif(temp.lt.40)then     
81        if (strictboundcia) then
82          print*,'Your temperatures are too low for this H2-H2 CIA dataset. If you '
83          print*,'really want to run simulations with hydrogen at T < 40 K, contact'
84          print*,'Robin Wordsworth [rwordsworth@uchicago.edu].'
85          stop
86        else
87          print*,'Your temperatures are too low for this H2-H2 CIA dataset'
88          print*,'you have chosen strictboundcia = ', strictboundcia
89          print*,'*********************************************************'
90          print*,' we allow model to continue but with temp = 40           '
91          print*,'  ...       for H2-H2 CIA dataset          ...           '
92          print*,'  ... we assume we know what you are doing ...           '
93          print*,'*********************************************************'
94          temp = 40
95        endif           
96      endif
97
98      amagat = (273.15/temp)*(pres/101325.0)
99
100      if(firstcall)then ! called by sugas_corrk only
101         print*,'----------------------------------------------------'
102         print*,'Initialising H2-H2 continuum from HITRAN database...'
103
104!     1.1 Open the ASCII files and set nS according to version
105         ! Only two possible versions for now : 2011 or 2018 (sanity check in inifis_mod)
106         if (versH2H2cia.eq.2011) then
107           nS = 2428
108           if (H2orthopara_mixture.eq."normal") then
109             dt_file=TRIM(datadir)//'/continuum_data/H2-H2_norm_2011.cia'
110           else if (H2orthopara_mixture.eq."equilibrium") then
111             dt_file=TRIM(datadir)//'/continuum_data/H2-H2_eq_2011.cia'
112           endif
113         else if (versH2H2cia.eq.2018) then
114           nS = 9600
115           if (H2orthopara_mixture.eq."normal") then
116             dt_file=TRIM(datadir)//'/continuum_data/H2-H2_norm_2018.cia'
117           else if (H2orthopara_mixture.eq."equilibrium") then
118             dt_file=TRIM(datadir)//'/continuum_data/H2-H2_eq_2018.cia'
119           endif
120         endif
121
122         if(.not.allocated(wn_arr))  allocate(wn_arr(nS))
123         if(.not.allocated(abs_arr)) allocate(abs_arr(nS,nT))
124
125!$OMP MASTER
126         open(33,file=dt_file,form='formatted',status='old',iostat=ios)
127         if (ios.ne.0) then        ! file not found
128           write(*,*) 'Error from interpolateH2H2'
129           write(*,*) 'Data file ',trim(dt_file),' not found.'
130           write(*,*) 'Check that your path to datagcm:',trim(datadir)
131           write(*,*) 'is correct. You can change it in callphys.def with:'
132           write(*,*) 'datadir = /absolute/path/to/datagcm'
133           write(*,*) 'Also check that the continuum data is there.'
134           call abort
135         else
136         
137         if(versH2H2cia.eq.2011) then
138           write(*,*) '... You are using H2-H2 CIA from 2011 but you should use more recent data available on HITRAN !'
139           write(*,*) '... (Especially if you are running a giant planet atmosphere)'
140           write(*,*) '... Just find out the H2-H2 CIA from 2018, put it in your datadir and have a look at interpolateH2H2.F90 ! .'
141         endif
142
143            do iT=1,nT
144               
145               ! Only two possibles values for now : 2011 or 2018 (sanity check in inifis_mod)
146               if(versH2H2cia.eq.2011) then
147                 read(33,fmat11) bleh,blah,blah,nres,Ttemp
148               else if (versH2H2cia.eq.2018) then
149                 read(33,fmat18) bleh,bleh,bleh,blah,blah,nres,Ttemp
150               endif
151
152               if(nS.ne.nres)then
153                  print*,'Resolution given in file: ',trim(dt_file)
154                  print*,'is ',nres,', which does not match nS.'
155                  print*,'Please adjust nS value in interpolateH2H2.F90'
156                  stop
157               endif
158               temp_arr(iT)=Ttemp
159
160               do k=1,nS
161                  read(33,*) wn_arr(k),abs_arr(k,it)
162               end do
163
164            end do
165     
166         endif
167         close(33)
168!$OMP END MASTER
169!$OMP BARRIER
170
171         print*,'interpolateH2H2: At wavenumber ',wn,' cm^-1'
172         print*,'   temperature ',temp,' K'
173         print*,'   pressure ',pres,' Pa'
174
175      endif
176     
177           
178           
179         call bilinearbig(nS,nT,wn_arr,temp_arr,abs_arr,wn,temp,abcoef,ind)
180
181         !print*,'the absorption is ',abcoef,' cm^5 molecule^-2'
182         !print*,'or ',abcoef*losch**2,' cm^-1 amagat^-2'
183
184         abcoef=abcoef*losch**2*100.0*amagat**2 ! convert to m^-1
185
186         !print*,'We have ',amagat,' amagats of H2'
187         !print*,'So the absorption is ',abcoef,' m^-1'
188
189         ! unlike for Rayleigh scattering, we do not currently weight by the BB function
190         ! however our bands are normally thin, so this is no big deal.
191
192      return
193    end subroutine interpolateH2H2
Note: See TracBrowser for help on using the repository browser.