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

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

Using of is_master keyword to allow nly the master processor to write the messages and avoid too heavy output files

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