source: trunk/LMDZ.GENERIC/libf/phystd/interpolateN2N2.F90

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

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

File size: 4.7 KB
RevLine 
[878]1subroutine interpolateN2N2(wn,temp,pres,abcoef,firstcall,ind)
[716]2
3  !==================================================================
4  !     
5  !     Purpose
6  !     -------
7  !     Calculates the N2-N2 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 pres               ! pressure               (Pascals)
[878]26  integer :: ind
[716]27
28  ! output
29  double precision abcoef             ! absorption coefficient (m^-1)
30
31  integer nS,nT
32  parameter(nS=582)
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 amagat
41  double precision wn_arr(nS)
42  double precision temp_arr(nT)
43  double precision abs_arr(nS,nT)
44
45  integer k,iT
46  logical firstcall
47
[1315]48  save wn_arr, temp_arr, abs_arr !read by master
[716]49
50  character*100 dt_file
51  integer strlen,ios
52
53  character(LEN=*), parameter :: fmat1 = "(A20,F10.3,F10.3,I7,F7.1,E10.3,F5.3)"
54
55  character*20 bleh
[2870]56  double precision blah, Ttemp, ztemp
[716]57  integer nres
58
[2870]59  ztemp=temp
[716]60
[2870]61  if(ztemp.gt.400)then
62   if (strictboundcia) then
[2662]63     if (is_master) then
[2870]64       print*,'Your temperatures are too high for this N2-N2 CIA dataset. If you '
65       print*,'really want to run simulations with hydrogen at T > 400 K,'
66       print*,'find relevant data.'
67     endif !is_master
[716]68     stop
[2870]69   else
70     !if (is_master) then
71     !  print*,'Your temperatures are too high for this N2-?2 CIA dataset'
72     !  print*,'you have chosen strictboundcia = ', strictboundcia
73     !  print*,'*********************************************************'
74     !  print*,' we allow model to continue but with temp = 4000          '
75     !  print*,'  ...       for N2-N2 CIA dataset          ...           '
76     !  print*,'  ... we assume we know what you are doing ...           '
77     !  print*,'*********************************************************'
78     !endif !is_master
79     ztemp = 400.
80   endif !strictboundcia
[716]81  endif
82
[2870]83
[716]84  amagat = (273.15/temp)*(pres/101325.0)
85
86  if(firstcall)then ! called by sugas_corrk only
[2662]87     if (is_master) print*,'----------------------------------------------------'
88     if (is_master) print*,'Initialising N2-N2 continuum from HITRAN database...'
[716]89
90     !     1.1 Open the ASCII files
91     dt_file=TRIM(datadir)//'/continuum_data/N2-N2_2011.cia'
92
[1315]93!$OMP MASTER
[716]94     open(33,file=dt_file,form='formatted',status='old',iostat=ios)
95     if (ios.ne.0) then        ! file not found
[2662]96        if (is_master) then
97          write(*,*) 'Error from interpolateN2N2'
98          write(*,*) 'Data file ',trim(dt_file),' not found.'
99          write(*,*) 'Check that your path to datagcm:',trim(datadir)
100          write(*,*) 'is correct. You can change it in callphys.def with:'
101          write(*,*) 'datadir = /absolute/path/to/datagcm'
102          write(*,*) 'Also check that the continuum data continuum_data/N2-N2_norm_2011.cia is there.'
103        endif
[716]104        call abort
105     else
106
107        do iT=1,nT
108
109           read(33,fmat1) bleh,blah,blah,nres,Ttemp
110           if(nS.ne.nres)then
[2662]111              if (is_master) then
112                print*,'Resolution given in file: ',trim(dt_file)
113                print*,'is ',nres,', which does not match nS.'
114                print*,'Please adjust nS value in interpolateN2N2.F90'
115              endif
[716]116              stop
117           endif
118           temp_arr(iT)=Ttemp
119
120           do k=1,nS
121              read(33,*) wn_arr(k),abs_arr(k,it)
122           end do
123
124        end do
125
126     endif
127     close(33)
[1315]128!$OMP END MASTER
129!$OMP BARRIER
[2662]130     if (is_master) then
131       print*,'interpolateN2N2: At wavenumber ',wn,' cm^-1'
132       print*,'   temperature ',temp,' K'
133       print*,'   pressure ',pres,' Pa'
134     endif
[878]135  endif
[2870]136     call bilinearbig(nS,nT,wn_arr,temp_arr,abs_arr,wn,ztemp,abcoef,ind)
[716]137
[878]138!     print*,'the absorption is ',abcoef,' cm^5 molecule^-2'
139!     print*,'or ',abcoef*losch**2,' cm^-1 amagat^-2'
[716]140
141     abcoef=abcoef*losch**2*100.0*amagat**2 ! convert to m^-1
[878]142!     abcoef=0.
[716]143
[878]144!     print*,'We have ',amagat,' amagats of N2'
145!     print*,'So the absorption is ',abcoef,' m^-1'
[716]146
[878]147!     unlike for Rayleigh scattering, we do not currently weight by the BB function
148!     however our bands are normally thin, so this is no big deal.
[716]149
150
151  return
152end subroutine interpolateN2N2
153
Note: See TracBrowser for help on using the repository browser.