source: trunk/LMDZ.GENERIC/libf/phystd/interpolateN2H2.F90 @ 3529

Last change on this file since 3529 was 2870, checked in by jleconte, 23 months ago

Removed too many prints in CIA interpolation + correction for hot H2H2

File size: 5.0 KB
RevLine 
[879]1subroutine interpolateN2H2(wn,temp,presN2,presH2,abcoef,firstcall,ind)
[716]2
3  !==================================================================
4  !     
5  !     Purpose
6  !     -------
7  !     Calculates the N2-H2 CIA opacity, using a lookup table from
8  !     HITRAN (2011)
9  !
10  !     Authors
11  !     -------
12  !     R. Wordsworth (2011)
13  !     
14  !==================================================================
15
[2870]16  use callkeys_mod, only: strictboundcia
[716]17  use datafile_mod, only: datadir
[2662]18  use mod_phys_lmdz_para, only : is_master
19
[716]20  implicit none
21
22  ! input
23  double precision wn                 ! wavenumber             (cm^-1)
24  double precision temp               ! temperature            (Kelvin)
25  double precision presN2             ! N2 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=1914)
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 amagatN2
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
[1315]49  save wn_arr, temp_arr, abs_arr !read by master
[716]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
[2870]57  double precision blah, Ttemp, ztemp
[716]58  integer nres
[879]59  integer ind
[716]60
[2870]61  ztemp=temp
62
63  if(ztemp.gt.400)then
64    if (strictboundcia) then
65      if (is_master) then
66        print*,'Your temperatures are too high for this N2-H2 CIA dataset. If you '
67        print*,'really want to run simulations with hydrogen at T > 400 K,'
68        print*,'find relevant data.'
69      endif !is_master
70      stop
71    else
72      !if (is_master) then
73      !  print*,'Your temperatures are too high for this N2-H2 CIA dataset'
74      !  print*,'you have chosen strictboundcia = ', strictboundcia
75      !  print*,'*********************************************************'
76      !  print*,' we allow model to continue but with temp = 4000          '
77      !  print*,'  ...       for N2-H2 CIA dataset          ...           '
78      !  print*,'  ... we assume we know what you are doing ...           '
79      !  print*,'*********************************************************'
80      !endif !is_master
81      ztemp = 400.
82    endif !strictboundcia
[716]83  endif
84
85  amagatN2 = (273.15/temp)*(presN2/101325.0)
86  amagatH2 = (273.15/temp)*(presH2/101325.0)
87
88  if(firstcall)then ! called by sugas_corrk only
[2662]89     if (is_master) print*,'----------------------------------------------------'
90     if (is_master) print*,'Initialising N2-H2 continuum from HITRAN database...'
[716]91
92     !     1.1 Open the ASCII files
93     dt_file=TRIM(datadir)//'/continuum_data/N2-H2_2011.cia'
94
[1315]95!$OMP MASTER
[716]96     open(33,file=dt_file,form='formatted',status='old',iostat=ios)
97     if (ios.ne.0) then        ! file not found
[2662]98        if (is_master) then
99          write(*,*) 'Error from interpolateN2H2'
100          write(*,*) 'Data file ',trim(dt_file),' not found.'
101          write(*,*) 'Check that your path to datagcm:',trim(datadir)
102          write(*,*) 'is correct. You can change it in callphys.def with:'
103          write(*,*) 'datadir = /absolute/path/to/datagcm'
104          write(*,*) 'Also check that the continuum data continuum_data/N2-H2_2011.cia is there.'
105        endif
[716]106        call abort
107     else
108
109        do iT=1,nT
110
111           read(33,fmat1) bleh,blah,blah,nres,Ttemp
112           if(nS.ne.nres)then
[2662]113              if (is_master) then
114                print*,'Resolution given in file: ',trim(dt_file)
115                print*,'is ',nres,', which does not match nS.'
116                print*,'Please adjust nS value in interpolateN2H2.F90'
117              endif
[716]118              stop
119           endif
120           temp_arr(iT)=Ttemp
121
122           do k=1,nS
123              read(33,*) wn_arr(k),abs_arr(k,it)
124           end do
125
126        end do
127
128     endif
129     close(33)
[1315]130!$OMP END MASTER
131!$OMP BARRIER
[2662]132     if (is_master) then
133       print*,'interpolateN2H2: At wavenumber ',wn,' cm^-1'
134       print*,'   temperature                 ',temp,' K'
135       print*,'   N2 partial pressure         ',presN2,' Pa'
136       print*,'   and H2 partial pressure     ',presH2,' Pa'
137     endif
[879]138  endif
[716]139
[2870]140  call bilinearbig(nS,nT,wn_arr,temp_arr,abs_arr,wn,ztemp,abcoef,ind)
[716]141
[879]142!     print*,'the absorption is ',abcoef,' cm^5 molecule^-2'
143!     print*,'or ',abcoef*losch**2,' cm^-1 amagat^-2'
[716]144
[879]145  abcoef=abcoef*losch**2*100.0*amagatN2*amagatH2 ! convert to m^-1
[716]146
[879]147!     print*,'We have ',amagatN2,' amagats of N2'
148!     print*,'and     ',amagatH2,' amagats of H2'
149!     print*,'So the absorption is ',abcoef,' m^-1'
[716]150
151
[879]152!  unlike for Rayleigh scattering, we do not currently weight by the BB function
153!  however our bands are normally thin, so this is no big deal.
[716]154
155
156  return
157end subroutine interpolateN2H2
158
Note: See TracBrowser for help on using the repository browser.