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

Last change on this file since 2719 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: 3.9 KB
Line 
1subroutine interpolateN2N2(wn,temp,pres,abcoef,firstcall,ind)
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
16  use datafile_mod, only: datadir
17  use mod_phys_lmdz_para, only : is_master
18
19  implicit none
20
21  ! input
22  double precision wn                 ! wavenumber             (cm^-1)
23  double precision temp               ! temperature            (Kelvin)
24  double precision pres               ! pressure               (Pascals)
25  integer :: ind
26
27  ! output
28  double precision abcoef             ! absorption coefficient (m^-1)
29
30  integer nS,nT
31  parameter(nS=582)
32  parameter(nT=10)
33
34  double precision, parameter :: losch = 2.6867774e19
35  ! Loschmit's number (molecule cm^-3 at STP)
36  ! converts cm^5 molecule^-2 --> cm^-1 amagat^-2
37  ! see Richard et al. 2011, JQSRT for details
38
39  double precision amagat
40  double precision wn_arr(nS)
41  double precision temp_arr(nT)
42  double precision abs_arr(nS,nT)
43
44  integer k,iT
45  logical firstcall
46
47  save wn_arr, temp_arr, abs_arr !read by master
48
49  character*100 dt_file
50  integer strlen,ios
51
52  character(LEN=*), parameter :: fmat1 = "(A20,F10.3,F10.3,I7,F7.1,E10.3,F5.3)"
53
54  character*20 bleh
55  double precision blah, Ttemp
56  integer nres
57
58
59  if(temp.gt.400)then
60     if (is_master) then
61       print*,'Your temperatures are too high for this N2-N2 CIA dataset.'
62       print*,'Currently, HITRAN provides data for this pair in the range 40-400 K.'
63     endif
64     stop
65  endif
66
67  amagat = (273.15/temp)*(pres/101325.0)
68
69  if(firstcall)then ! called by sugas_corrk only
70     if (is_master) print*,'----------------------------------------------------'
71     if (is_master) print*,'Initialising N2-N2 continuum from HITRAN database...'
72
73     !     1.1 Open the ASCII files
74     dt_file=TRIM(datadir)//'/continuum_data/N2-N2_2011.cia'
75
76!$OMP MASTER
77     open(33,file=dt_file,form='formatted',status='old',iostat=ios)
78     if (ios.ne.0) then        ! file not found
79        if (is_master) then
80          write(*,*) 'Error from interpolateN2N2'
81          write(*,*) 'Data file ',trim(dt_file),' not found.'
82          write(*,*) 'Check that your path to datagcm:',trim(datadir)
83          write(*,*) 'is correct. You can change it in callphys.def with:'
84          write(*,*) 'datadir = /absolute/path/to/datagcm'
85          write(*,*) 'Also check that the continuum data continuum_data/N2-N2_norm_2011.cia is there.'
86        endif
87        call abort
88     else
89
90        do iT=1,nT
91
92           read(33,fmat1) bleh,blah,blah,nres,Ttemp
93           if(nS.ne.nres)then
94              if (is_master) then
95                print*,'Resolution given in file: ',trim(dt_file)
96                print*,'is ',nres,', which does not match nS.'
97                print*,'Please adjust nS value in interpolateN2N2.F90'
98              endif
99              stop
100           endif
101           temp_arr(iT)=Ttemp
102
103           do k=1,nS
104              read(33,*) wn_arr(k),abs_arr(k,it)
105           end do
106
107        end do
108
109     endif
110     close(33)
111!$OMP END MASTER
112!$OMP BARRIER
113     if (is_master) then
114       print*,'interpolateN2N2: At wavenumber ',wn,' cm^-1'
115       print*,'   temperature ',temp,' K'
116       print*,'   pressure ',pres,' Pa'
117     endif
118  endif
119     call bilinearbig(nS,nT,wn_arr,temp_arr,abs_arr,wn,temp,abcoef,ind)
120
121!     print*,'the absorption is ',abcoef,' cm^5 molecule^-2'
122!     print*,'or ',abcoef*losch**2,' cm^-1 amagat^-2'
123
124     abcoef=abcoef*losch**2*100.0*amagat**2 ! convert to m^-1
125!     abcoef=0.
126
127!     print*,'We have ',amagat,' amagats of N2'
128!     print*,'So the absorption is ',abcoef,' m^-1'
129
130!     unlike for Rayleigh scattering, we do not currently weight by the BB function
131!     however our bands are normally thin, so this is no big deal.
132
133
134  return
135end subroutine interpolateN2N2
136
Note: See TracBrowser for help on using the repository browser.