source: trunk/LMDZ.GENERIC/libf/phystd/interpolateH2CH4.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

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