source: trunk/LMDZ.GENERIC/libf/phystd/interpolateHeCH4.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.

File size: 5.3 KB
Line 
1subroutine interpolateHeCH4(wn,temp,presHe,presCH4,abcoef,firstcall,ind)
2
3  !==================================================================
4  !     
5  !     Purpose
6  !     -------
7  !     Calculates the He-CH4 CIA opacity, using a lookup table from
8  !     HITRAN (2018)
9  !
10  !     Authors
11  !     -------
12  !     R. Wordsworth (2011)
13  !     
14  !==================================================================
15
16  use callkeys_mod, only: strictboundcia
17  use datafile_mod, only: datadir
18  implicit none
19
20  ! input
21  double precision wn                 ! wavenumber             (cm^-1)
22  double precision temp               ! temperature            (Kelvin)
23  double precision presCH4            ! CH4 partial pressure    (Pascals)
24  double precision presHe             ! He partial pressure    (Pascals)
25
26  ! output
27  double precision abcoef             ! absorption coefficient (m^-1)
28
29  integer nS,nT
30  parameter(nS=1000)
31  parameter(nT=10)
32
33  double precision, parameter :: losch = 2.6867774e19
34  ! Loschmit's number (molecule cm^-3 at STP)
35  ! converts cm^5 molecule^-2 --> cm^-1 amagat^-2
36  ! see Richard et al. 2011, JQSRT for details
37
38  double precision amagatCH4
39  double precision amagatHe
40  double precision wn_arr(nS)
41  double precision temp_arr(nT)
42  double precision abs_arr(nS,nT)
43
44  integer k,iT
45  logical firstcall
46
47  save wn_arr, temp_arr, abs_arr !read by master
48
49  character*100 dt_file
50  integer strlen,ios
51
52  character(LEN=*), parameter :: fmat1 = "(A20,F10.3,F10.3,I7,F7.1,E10.3,F5.3)"
53
54  character*20 bleh
55  double precision blah, Ttemp
56  integer nres
57  integer ind
58
59  if(temp.gt.350)then
60    if (strictboundcia) then
61      print*,'Your temperatures are too high for this He-CH4 CIA dataset.'
62      print*,'Please run mixed He-CH4 atmospheres below T = 350 K.'     
63      stop
64    else
65      print*,'Your temperatures are too high for this He-CH4 CIA dataset'
66      print*,'you have chosen strictboundcia = ', strictboundcia
67      print*,'*********************************************************'
68      print*,' we allow model to continue but with temp = 350          '
69      print*,'  ...       for He-CH4 CIA dataset         ...           '
70      print*,'  ... we assume we know what you are doing ...           '
71      print*,'*********************************************************'
72      temp = 350
73    endif
74  elseif(temp.lt.40)then
75    if (strictboundcia) then
76      print*,'Your temperatures are too low for this He-CH4 CIA dataset.'
77      print*,'Please run mixed He-CH4 atmospheres above T = 40 K.'     
78      stop
79    else
80      print*,'Your temperatures are too low for this He-CH4 CIA dataset'
81      print*,'you have chosen strictboundcia = ', strictboundcia
82      print*,'*********************************************************'
83      print*,' we allow model to continue but with temp = 40           '
84      print*,'  ...       for He-CH4 CIA dataset         ...           '
85      print*,'  ... we assume we know what you are doing ...           '
86      print*,'*********************************************************'
87      temp = 40
88    endif
89  endif
90
91  amagatCH4 = (273.15/temp)*(presCH4/101325.0)
92  amagatHe = (273.15/temp)*(presHe/101325.0)
93
94  if(firstcall)then ! called by sugas_corrk only
95     print*,'----------------------------------------------------'
96     print*,'Initialising He-CH4 continuum from HITRAN database...'
97
98     !     1.1 Open the ASCII files
99     dt_file=TRIM(datadir)//'/continuum_data/CH4-He_2018.cia'
100
101!$OMP MASTER
102     open(33,file=dt_file,form='formatted',status='old',iostat=ios)
103     if (ios.ne.0) then        ! file not found
104        write(*,*) 'Error from interpolateHeCH4'
105        write(*,*) 'Data file ',trim(dt_file),' not found.'
106        write(*,*) 'Check that your path to datagcm:',trim(datadir)
107        write(*,*) 'is correct. You can change it in callphys.def with:'
108        write(*,*) 'datadir = /absolute/path/to/datagcm'
109        write(*,*) 'Also check that the continuum data continuum_data/He-CH4_2018.cia is there.'
110        call abort
111     else
112
113        do iT=1,nT
114
115           read(33,fmat1) bleh,blah,blah,nres,Ttemp
116           if(nS.ne.nres)then
117              print*,'Resolution given in file: ',trim(dt_file)
118              print*,'is ',nres,', which does not match nS.'
119              print*,'Please adjust nS value in interpolateHeCH4.F90'
120              stop
121           endif
122           temp_arr(iT)=Ttemp
123
124           do k=1,nS
125              read(33,*) wn_arr(k),abs_arr(k,it)
126           end do
127
128        end do
129
130     endif
131     close(33)
132!$OMP END MASTER
133!$OMP BARRIER
134
135     print*,'interpolateHeCH4: At wavenumber ',wn,' cm^-1'
136     print*,'   temperature                 ',temp,' K'
137     print*,'   CH4 partial pressure        ',presCH4,' Pa'
138     print*,'   and He partial pressure     ',presHe,' Pa'
139
140  endif
141
142  call bilinearbig(nS,nT,wn_arr,temp_arr,abs_arr,wn,temp,abcoef,ind)
143
144!     print*,'the absorption is ',abcoef,' cm^5 molecule^-2'
145!     print*,'or ',abcoef*losch**2,' cm^-1 amagat^-2'
146
147  abcoef=abcoef*losch**2*100.0*amagatCH4*amagatHe ! convert to m^-1
148
149!     print*,'We have ',amagatCH4,' amagats of CH4'
150!     print*,'and     ',amagatHe,' amagats of He'
151!     print*,'So the absorption is ',abcoef,' m^-1'
152
153
154!  unlike for Rayleigh scattering, we do not currently weight by the BB function
155!  however our bands are normally thin, so this is no big deal.
156
157
158  return
159end subroutine interpolateHeCH4
160
Note: See TracBrowser for help on using the repository browser.