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

Last change on this file since 848 was 789, checked in by aslmd, 12 years ago

LMDZ.GENERIC
A more robust way to count lines in setspi and setspv.
bandlen.txt file is no longer used. This was causing problems with MPI computations.

File size: 6.9 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
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=-1  !because first line is not an actual value
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=.false.
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      do while (ierr==0)
97        read(131,*,iostat=ierr) file_entries
98        if (ierr==0) nb=nb+1
99      enddo
100      close(131)
101      write(*,*) 'setspi: L_NSPECTI = ',L_NSPECTI, 'in the model '
102      write(*,*) '        there are   ',nb, 'entries in ',TRIM(file_path)
103      if(nb.ne.L_NSPECTI) then
104         write(*,*) 'MISMATCH !! I stop here'
105         call abort
106      endif
107
108      ! load and display the data
109      open(111,file=file_path,form='formatted')
110      read(111,*)
111      do M=1,L_NSPECTI-1
112         read(111,*) BWNI(M)
113      end do
114      read(111,*) lastband
115      close(111)
116      BWNI(L_NSPECTI)  =lastband(1)
117      BWNI(L_NSPECTI+1)=lastband(2)
118
119      print*,''
120      print*,'setspi: IR band limits:'
121      do M=1,L_NSPECTI+1
122         print*,m,'-->',BWNI(M),' cm^-1'
123      end do
124
125!     Set up mean wavenumbers and wavenumber deltas.  Units of
126!     wavenumbers is cm^(-1); units of wavelengths is microns.
127
128      do M=1,L_NSPECTI
129         WNOI(M)  = 0.5*(BWNI(M+1)+BWNI(M))
130         DWNI(M)  = BWNI(M+1)-BWNI(M)
131         WAVEI(M) = 1.0E+4/WNOI(M)
132         BLAMI(M) = 0.01/BWNI(M)         
133      end do
134      BLAMI(M) = 0.01/BWNI(M)
135!     note M=L_NSPECTI+1 after loop due to Fortran bizarreness
136
137!=======================================================================
138!     For each IR wavelength interval, compute the integral of B(T), the
139!     Planck function, divided by the wavelength interval, in cm-1.  The
140!     integration is in MKS units, the final answer is the same as the
141!     original planck.f; W m^-2 wavenumber^-1, where wavenumber is in CM^-1.
142
143      print*,''
144      print*,'setspi: Current Planck integration range:'
145      print*,'T = ',dble(NTstar)/NTfac, ' to ',dble(NTstop)/NTfac,' K.'
146
147      do NW=1,L_NSPECTI
148         a = 1.0D-2/BWNI(NW+1)
149         b = 1.0D-2/BWNI(NW)
150         bpa = (b+a)/2.0
151         bma = (b-a)/2.0
152         do nt=NTstar,NTstop
153            T   = dble(NT)/NTfac
154            ans = 0.0D0
155
156            do mm=1,12
157               y    = bma*x(mm)+bpa
158               ans  = ans + w(mm)*c1/(y**5*(exp(c2/(y*T))-1.0D0))
159            end do
160
161            planckir(NW,nt-NTstar+1) = ans*bma/(PI*DWNI(NW))
162         end do
163      end do
164         
165      ! force planck=sigma*eps*T^4 for each temperature in array
166      if(forceEC)then
167         print*,'setspi: Force F=sigma*eps*T^4 for all values of T!'
168         do nt=NTstar,NTstop
169            plancksum=0.0
170            T=dble(NT)/NTfac
171       
172            do NW=1,L_NSPECTI
173               plancksum=plancksum+  &
174                  planckir(NW,nt-NTstar+1)*DWNI(NW)*pi
175            end do
176
177            do NW=1,L_NSPECTI
178               planckir(NW,nt-NTstar+1)=     &
179                  planckir(NW,nt-NTstar+1)*  &
180                          sigma*(dble(nt)/NTfac)**4/plancksum
181            end do
182         end do
183      endif
184
185      if(planckcheck)then
186         ! check energy conservation at lower temperature boundary
187         plancksum=0.0
188         nt=NTstar
189         do NW=1,L_NSPECTI
190            plancksum=plancksum+planckir(NW,nt-NTstar+1)*DWNI(NW)*pi
191         end do
192         print*,'setspi: At lower limit:'
193         print*,'in model sig*T^4 = ',plancksum,' W m^-2'
194         print*,'actual sig*T^4   = ',sigma*(dble(nt)/NTfac)**4,' W m^-2'
195         
196         ! check energy conservation at upper temperature boundary
197         plancksum=0.0
198         nt=NTstop
199         do NW=1,L_NSPECTI
200            plancksum=plancksum+planckir(NW,nt-NTstar+1)*DWNI(NW)*pi
201         end do
202         print*,'setspi: At upper limit:'
203         print*,'in model sig*T^4 = ',plancksum,' W m^-2'
204         print*,'actual sig*T^4   = ',sigma*(dble(nt)/NTfac)**4,' W m^-2'
205         print*,''
206      endif
207
208      return
209    end subroutine setspi
Note: See TracBrowser for help on using the repository browser.