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

Last change on this file since 2633 was 2633, checked in by dbardet, 3 years ago

Adding of strictboundH2H2cia flag to allow the run to continue even for temperature bigger than 400K (that is out of the current H2H2 cia tables)

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