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

Last change on this file since 2613 was 2245, checked in by jvatant, 5 years ago

+ Add a 'versH2H2cia' int key in callphys that allows two values (2011 or 2018) to

deal with updated HITRAN file (for interpolateH2H2.F90) from 2018 that includes the
H2H2 dimer from Fletcher et al. 2018, useful for giant planets.
Retrocompatibility is ok, default value to 2011.

--JVO

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