source: trunk/LMDZ.GENERIC/libf/phystd/setspi.F90 @ 989

Last change on this file since 989 was 989, checked in by emillour, 11 years ago

Generic GCM:

  • Some minor changes so that gcm compiles with gfortran:
    • Added option to compile "long lines" (>132 characters) in makegcm_gfortran
    • Removed use of isnan() in physiq.F90 (it is not a standard function)
    • Avoid possible underflow of psat in watercommon_h.F90
    • Adapted the checks on the *IR and *VI band files to be more strict

EM

File size: 7.0 KB
Line 
1      subroutine setspi
2
3!==================================================================
4!     
5!     Purpose
6!     -------
7!     Set up spectral intervals and Planck function in the longwave.
8!     
9!     Authors
10!     -------
11!     Adapted from setspi in the NASA Ames radiative code by
12!     Robin Wordsworth (2009).
13!     
14!     Called by
15!     ---------
16!     callcorrk.F
17!     
18!     Calls
19!     -----
20!     none
21!     
22!==================================================================
23
24      use radinc_h,    only: L_NSPECTI,corrkdir,banddir,NTstar,NTstop,NTfac
25      use radcommon_h, only: BWNI,BLAMI,WNOI,DWNI,WAVEI,planckir,sigma
26      use datafile_mod, only: datadir
27
28      implicit none
29
30#include "callkeys.h"
31#include "comcstfi.h"
32
33      logical file_ok
34      integer nw, nt, m, mm, file_entries
35      real*8 a, b, ans, y, bpa, bma, T, dummy
36
37      character(len=30)  :: temp1
38      character(len=200) :: file_id
39      character(len=200) :: file_path
40
41!     C1 and C2 values from Goody and Yung (2nd edition)  MKS units
42!     These values lead to a "sigma" (sigma*T^4) of 5.67032E-8 W m^-2 K^-4
43
44      real*8 :: c1 = 3.741832D-16 ! W m^-2
45      real*8 :: c2 = 1.438786D-2  ! m K
46     
47      real*8 :: lastband(2), plancksum
48
49      !! used to count lines
50      integer :: nb=0
51      integer :: ierr=0
52
53      logical forceEC, planckcheck
54
55      real*8 :: x(12) = [ -0.981560634246719D0,  -0.904117256370475D0, &
56      -0.769902674194305D0,  -0.587317954286617D0,                     &
57      -0.367831498998180D0,  -0.125233408511469D0,                     &
58       0.125233408511469D0,   0.367831498998180D0,                     &
59       0.587317954286617D0,   0.769902674194305D0,                     &
60       0.904117256370475D0,   0.981560634246719D0  ]
61
62      real*8 :: w(12) = [  0.047175336386512D0,   0.106939325995318D0, &
63           0.160078328543346D0,   0.203167426723066D0,                 &
64           0.233492536538355D0,   0.249147045813403D0,                 &
65           0.249147045813403D0,   0.233492536538355D0,                 &
66           0.203167426723066D0,   0.160078328543346D0,                 &
67           0.106939325995318D0,   0.047175336386512D0  ]
68      mm=0
69
70      forceEC=.true.
71      planckcheck=.true.
72
73!=======================================================================
74!     Set up spectral bands - wavenumber [cm^(-1)]. Go from smaller to
75!     larger wavenumbers.
76
77      write(temp1,'(i2.2)') L_NSPECTI
78      !file_id='/corrk_data/' // corrkdir(1:LEN_TRIM(corrkdir)) // '/narrowbands_IR.in'
79      file_id='/corrk_data/'//trim(adjustl(banddir))//'/narrowbands_IR.in'
80      file_path=TRIM(datadir)//TRIM(file_id)
81
82      ! check that the file exists
83      inquire(FILE=file_path,EXIST=file_ok)
84      if(.not.file_ok) then
85         write(*,*)'The file ',TRIM(file_path)
86         write(*,*)'was not found by setspi.F90, exiting.'
87         write(*,*)'Check that your path to datagcm:',trim(datadir)
88         write(*,*)' is correct. You can change it in callphys.def with:'
89         write(*,*)' datadir = /absolute/path/to/datagcm'
90         write(*,*)'Also check that the corrkdir you chose in callphys.def exists.'
91         call abort
92      endif
93   
94      ! check that the file contains the right number of bands
95      open(131,file=file_path,form='formatted')
96      read(131,*,iostat=ierr) file_entries
97      do while (ierr==0)
98        read(131,*,iostat=ierr) dummy
99!        write(*,*) 'setspi: file_entries:',dummy,'ierr=',ierr
100        if (ierr==0) nb=nb+1
101      enddo
102      close(131)
103
104      write(*,*) 'setspi: L_NSPECTI = ',L_NSPECTI, 'in the model '
105      write(*,*) '        there are   ',nb, 'entries in ',TRIM(file_path)
106      if(nb.ne.L_NSPECTI) then
107         write(*,*) 'MISMATCH !! I stop here'
108         call abort
109      endif
110
111      ! load and display the data
112      open(111,file=file_path,form='formatted')
113      read(111,*)
114      do M=1,L_NSPECTI-1
115         read(111,*) BWNI(M)
116      end do
117      read(111,*) lastband
118      close(111)
119      BWNI(L_NSPECTI)  =lastband(1)
120      BWNI(L_NSPECTI+1)=lastband(2)
121
122      print*,''
123      print*,'setspi: IR band limits:'
124      do M=1,L_NSPECTI+1
125         print*,m,'-->',BWNI(M),' cm^-1'
126      end do
127
128!     Set up mean wavenumbers and wavenumber deltas.  Units of
129!     wavenumbers is cm^(-1); units of wavelengths is microns.
130
131      do M=1,L_NSPECTI
132         WNOI(M)  = 0.5D0*(BWNI(M+1)+BWNI(M))
133         DWNI(M)  = BWNI(M+1)-BWNI(M)
134         WAVEI(M) = 1.0D+4/WNOI(M)
135         BLAMI(M) = 0.01D0/BWNI(M)         
136      end do
137      BLAMI(M) = 0.01D0/BWNI(M)
138!     note M=L_NSPECTI+1 after loop due to Fortran bizarreness
139
140!=======================================================================
141!     For each IR wavelength interval, compute the integral of B(T), the
142!     Planck function, divided by the wavelength interval, in cm-1.  The
143!     integration is in MKS units, the final answer is the same as the
144!     original planck.f; W m^-2 wavenumber^-1, where wavenumber is in CM^-1.
145
146      print*,''
147      print*,'setspi: Current Planck integration range:'
148      print*,'T = ',dble(NTstar)/NTfac, ' to ',dble(NTstop)/NTfac,' K.'
149
150      do NW=1,L_NSPECTI
151         a = 1.0D-2/BWNI(NW+1)
152         b = 1.0D-2/BWNI(NW)
153         bpa = (b+a)/2.0D0
154         bma = (b-a)/2.0D0
155         do nt=NTstar,NTstop
156            T   = dble(NT)/NTfac
157            ans = 0.0D0
158
159            do mm=1,12
160               y    = bma*x(mm)+bpa
161               ans  = ans + w(mm)*c1/(y**5*(exp(c2/(y*T))-1.0D0))
162            end do
163
164            planckir(NW,nt-NTstar+1) = ans*bma/(PI*DWNI(NW))
165         end do
166      end do
167         
168      ! force planck=sigma*eps*T^4 for each temperature in array
169      if(forceEC)then
170         print*,'setspi: Force F=sigma*eps*T^4 for all values of T!'
171         do nt=NTstar,NTstop
172            plancksum=0.0D0
173            T=dble(NT)/NTfac
174       
175            do NW=1,L_NSPECTI
176               plancksum=plancksum+  &
177                  planckir(NW,nt-NTstar+1)*DWNI(NW)*pi
178            end do
179
180            do NW=1,L_NSPECTI
181               planckir(NW,nt-NTstar+1)=     &
182                  planckir(NW,nt-NTstar+1)*  &
183                          sigma*(dble(nt)/NTfac)**4/plancksum
184            end do
185         end do
186      endif
187
188      if(planckcheck)then
189         ! check energy conservation at lower temperature boundary
190         plancksum=0.0D0
191         nt=NTstar
192         do NW=1,L_NSPECTI
193            plancksum=plancksum+planckir(NW,nt-NTstar+1)*DWNI(NW)*pi
194         end do
195         print*,'setspi: At lower limit:'
196         print*,'in model sig*T^4 = ',plancksum,' W m^-2'
197         print*,'actual sig*T^4   = ',sigma*(dble(nt)/NTfac)**4,' W m^-2'
198         
199         ! check energy conservation at upper temperature boundary
200         plancksum=0.0D0
201         nt=NTstop
202         do NW=1,L_NSPECTI
203            plancksum=plancksum+planckir(NW,nt-NTstar+1)*DWNI(NW)*pi
204         end do
205         print*,'setspi: At upper limit:'
206         print*,'in model sig*T^4 = ',plancksum,' W m^-2'
207         print*,'actual sig*T^4   = ',sigma*(dble(nt)/NTfac)**4,' W m^-2'
208         print*,''
209      endif
210
211      return
212    end subroutine setspi
Note: See TracBrowser for help on using the repository browser.