source: trunk/LMDZ.GENERIC/libf/phystd/interpolateH2CH4.F90 @ 2870

Last change on this file since 2870 was 2870, checked in by jleconte, 23 months ago

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

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