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

Last change on this file since 1243 was 997, checked in by emillour, 12 years ago

Generic GCM:

  • bug fix in setspi & setspv : counters must be initialized in routine in case of multiple calls

JL

File size: 7.0 KB
RevLine 
[135]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
[543]24      use radinc_h,    only: L_NSPECTI,corrkdir,banddir,NTstar,NTstop,NTfac
[135]25      use radcommon_h, only: BWNI,BLAMI,WNOI,DWNI,WAVEI,planckir,sigma
[374]26      use datafile_mod, only: datadir
[135]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
[989]35      real*8 a, b, ans, y, bpa, bma, T, dummy
[135]36
37      character(len=30)  :: temp1
[716]38      character(len=200) :: file_id
39      character(len=200) :: file_path
[135]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
[789]49      !! used to count lines
[997]50      integer :: nb
51      integer :: ierr
[789]52
[135]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
[959]70      forceEC=.true.
[135]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'
[374]80      file_path=TRIM(datadir)//TRIM(file_id)
[135]81
82      ! check that the file exists
83      inquire(FILE=file_path,EXIST=file_ok)
84      if(.not.file_ok) then
[374]85         write(*,*)'The file ',TRIM(file_path)
[135]86         write(*,*)'was not found by setspi.F90, exiting.'
[374]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.'
[135]91         call abort
92      endif
[789]93   
[997]94      nb=0
95      ierr=0
[789]96      ! check that the file contains the right number of bands
97      open(131,file=file_path,form='formatted')
[989]98      read(131,*,iostat=ierr) file_entries
[789]99      do while (ierr==0)
[989]100        read(131,*,iostat=ierr) dummy
101!        write(*,*) 'setspi: file_entries:',dummy,'ierr=',ierr
[789]102        if (ierr==0) nb=nb+1
103      enddo
[135]104      close(131)
[989]105
[789]106      write(*,*) 'setspi: L_NSPECTI = ',L_NSPECTI, 'in the model '
107      write(*,*) '        there are   ',nb, 'entries in ',TRIM(file_path)
108      if(nb.ne.L_NSPECTI) then
109         write(*,*) 'MISMATCH !! I stop here'
[135]110         call abort
111      endif
112
113      ! load and display the data
114      open(111,file=file_path,form='formatted')
115      read(111,*)
116      do M=1,L_NSPECTI-1
117         read(111,*) BWNI(M)
118      end do
119      read(111,*) lastband
120      close(111)
121      BWNI(L_NSPECTI)  =lastband(1)
122      BWNI(L_NSPECTI+1)=lastband(2)
123
124      print*,''
[374]125      print*,'setspi: IR band limits:'
[135]126      do M=1,L_NSPECTI+1
127         print*,m,'-->',BWNI(M),' cm^-1'
128      end do
129
130!     Set up mean wavenumbers and wavenumber deltas.  Units of
131!     wavenumbers is cm^(-1); units of wavelengths is microns.
132
133      do M=1,L_NSPECTI
[961]134         WNOI(M)  = 0.5D0*(BWNI(M+1)+BWNI(M))
[135]135         DWNI(M)  = BWNI(M+1)-BWNI(M)
[961]136         WAVEI(M) = 1.0D+4/WNOI(M)
137         BLAMI(M) = 0.01D0/BWNI(M)         
[135]138      end do
[961]139      BLAMI(M) = 0.01D0/BWNI(M)
[135]140!     note M=L_NSPECTI+1 after loop due to Fortran bizarreness
141
142!=======================================================================
143!     For each IR wavelength interval, compute the integral of B(T), the
144!     Planck function, divided by the wavelength interval, in cm-1.  The
145!     integration is in MKS units, the final answer is the same as the
146!     original planck.f; W m^-2 wavenumber^-1, where wavenumber is in CM^-1.
147
148      print*,''
[374]149      print*,'setspi: Current Planck integration range:'
[543]150      print*,'T = ',dble(NTstar)/NTfac, ' to ',dble(NTstop)/NTfac,' K.'
[135]151
152      do NW=1,L_NSPECTI
153         a = 1.0D-2/BWNI(NW+1)
154         b = 1.0D-2/BWNI(NW)
[961]155         bpa = (b+a)/2.0D0
156         bma = (b-a)/2.0D0
[135]157         do nt=NTstar,NTstop
[543]158            T   = dble(NT)/NTfac
[135]159            ans = 0.0D0
160
161            do mm=1,12
162               y    = bma*x(mm)+bpa
163               ans  = ans + w(mm)*c1/(y**5*(exp(c2/(y*T))-1.0D0))
164            end do
165
166            planckir(NW,nt-NTstar+1) = ans*bma/(PI*DWNI(NW))
167         end do
168      end do
169         
170      ! force planck=sigma*eps*T^4 for each temperature in array
171      if(forceEC)then
[374]172         print*,'setspi: Force F=sigma*eps*T^4 for all values of T!'
[135]173         do nt=NTstar,NTstop
[961]174            plancksum=0.0D0
[543]175            T=dble(NT)/NTfac
[135]176       
177            do NW=1,L_NSPECTI
178               plancksum=plancksum+  &
179                  planckir(NW,nt-NTstar+1)*DWNI(NW)*pi
180            end do
181
182            do NW=1,L_NSPECTI
183               planckir(NW,nt-NTstar+1)=     &
184                  planckir(NW,nt-NTstar+1)*  &
[543]185                          sigma*(dble(nt)/NTfac)**4/plancksum
[135]186            end do
187         end do
188      endif
189
190      if(planckcheck)then
191         ! check energy conservation at lower temperature boundary
[961]192         plancksum=0.0D0
[135]193         nt=NTstar
194         do NW=1,L_NSPECTI
195            plancksum=plancksum+planckir(NW,nt-NTstar+1)*DWNI(NW)*pi
196         end do
[374]197         print*,'setspi: At lower limit:'
[135]198         print*,'in model sig*T^4 = ',plancksum,' W m^-2'
[543]199         print*,'actual sig*T^4   = ',sigma*(dble(nt)/NTfac)**4,' W m^-2'
[135]200         
201         ! check energy conservation at upper temperature boundary
[961]202         plancksum=0.0D0
[135]203         nt=NTstop
204         do NW=1,L_NSPECTI
205            plancksum=plancksum+planckir(NW,nt-NTstar+1)*DWNI(NW)*pi
206         end do
[374]207         print*,'setspi: At upper limit:'
[135]208         print*,'in model sig*T^4 = ',plancksum,' W m^-2'
[543]209         print*,'actual sig*T^4   = ',sigma*(dble(nt)/NTfac)**4,' W m^-2'
[135]210         print*,''
211      endif
212
213      return
214    end subroutine setspi
Note: See TracBrowser for help on using the repository browser.