source: trunk/LMDZ.GENERIC/libf/phystd/interpolateHeCH4.F90 @ 2757

Last change on this file since 2757 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

File size: 5.7 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  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 presHe             ! He partial pressure    (Pascals)
27
28  ! output
29  double precision abcoef             ! absorption coefficient (m^-1)
30
31  integer nS,nT
32  parameter(nS=1000)
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 amagatHe
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
58  integer nres
59  integer ind
60
61  if(temp.gt.350)then
62    if (strictboundcia) then
63      if (is_master) then
64        print*,'Your temperatures are too high for this He-CH4 CIA dataset.'
65        print*,'Please run mixed He-CH4 atmospheres below T = 350 K.'     
66      endif
67      stop
68    else
69      if (is_master) then
70        print*,'Your temperatures are too high for this He-CH4 CIA dataset'
71        print*,'you have chosen strictboundcia = ', strictboundcia
72        print*,'*********************************************************'
73        print*,' we allow model to continue but with temp = 350          '
74        print*,'  ...       for He-CH4 CIA dataset         ...           '
75        print*,'  ... we assume we know what you are doing ...           '
76        print*,'*********************************************************'
77      endif
78      temp = 350
79    endif
80  elseif(temp.lt.40)then
81    if (strictboundcia) then
82      if (is_master) then
83        print*,'Your temperatures are too low for this He-CH4 CIA dataset.'
84        print*,'Please run mixed He-CH4 atmospheres above T = 40 K.'     
85      endif
86      stop
87    else
88      if (is_master) then
89        print*,'Your temperatures are too low for this He-CH4 CIA dataset'
90        print*,'you have chosen strictboundcia = ', strictboundcia
91        print*,'*********************************************************'
92        print*,' we allow model to continue but with temp = 40           '
93        print*,'  ...       for He-CH4 CIA dataset         ...           '
94        print*,'  ... we assume we know what you are doing ...           '
95        print*,'*********************************************************'
96      endif
97      temp = 40
98    endif
99  endif
100
101  amagatCH4 = (273.15/temp)*(presCH4/101325.0)
102  amagatHe = (273.15/temp)*(presHe/101325.0)
103
104  if(firstcall)then ! called by sugas_corrk only
105     if (is_master) print*,'----------------------------------------------------'
106     if (is_master) print*,'Initialising He-CH4 continuum from HITRAN database...'
107
108     !     1.1 Open the ASCII files
109     dt_file=TRIM(datadir)//'/continuum_data/CH4-He_2018.cia'
110
111!$OMP MASTER
112     open(33,file=dt_file,form='formatted',status='old',iostat=ios)
113     if (ios.ne.0) then        ! file not found
114        if (is_master) then
115          write(*,*) 'Error from interpolateHeCH4'
116          write(*,*) 'Data file ',trim(dt_file),' not found.'
117          write(*,*) 'Check that your path to datagcm:',trim(datadir)
118          write(*,*) 'is correct. You can change it in callphys.def with:'
119          write(*,*) 'datadir = /absolute/path/to/datagcm'
120          write(*,*) 'Also check that the continuum data continuum_data/He-CH4_2018.cia is there.'
121        endif
122        call abort
123     else
124
125        do iT=1,nT
126
127           read(33,fmat1) bleh,blah,blah,nres,Ttemp
128           if(nS.ne.nres)then
129              if (is_master) then
130                print*,'Resolution given in file: ',trim(dt_file)
131                print*,'is ',nres,', which does not match nS.'
132                print*,'Please adjust nS value in interpolateHeCH4.F90'
133              endif
134              stop
135           endif
136           temp_arr(iT)=Ttemp
137
138           do k=1,nS
139              read(33,*) wn_arr(k),abs_arr(k,it)
140           end do
141
142        end do
143
144     endif
145     close(33)
146!$OMP END MASTER
147!$OMP BARRIER
148     if (is_master) then
149       print*,'interpolateHeCH4: At wavenumber ',wn,' cm^-1'
150       print*,'   temperature                 ',temp,' K'
151       print*,'   CH4 partial pressure        ',presCH4,' Pa'
152       print*,'   and He partial pressure     ',presHe,' Pa'
153     endif
154  endif
155
156  call bilinearbig(nS,nT,wn_arr,temp_arr,abs_arr,wn,temp,abcoef,ind)
157
158!     print*,'the absorption is ',abcoef,' cm^5 molecule^-2'
159!     print*,'or ',abcoef*losch**2,' cm^-1 amagat^-2'
160
161  abcoef=abcoef*losch**2*100.0*amagatCH4*amagatHe ! convert to m^-1
162
163!     print*,'We have ',amagatCH4,' amagats of CH4'
164!     print*,'and     ',amagatHe,' amagats of He'
165!     print*,'So the absorption is ',abcoef,' m^-1'
166
167
168!  unlike for Rayleigh scattering, we do not currently weight by the BB function
169!  however our bands are normally thin, so this is no big deal.
170
171
172  return
173end subroutine interpolateHeCH4
174
Note: See TracBrowser for help on using the repository browser.