source: codes/icosagcm/branches/SATURN_DYNAMICO/LMDZ.COMMON/libf/phystd/setspi.F90

Last change on this file was 313, checked in by ymipsl, 10 years ago
  • implement splitting of XIOS file for lmdz physics
  • Termination is done properly in parallel by calling MPI_ABORT instead of abort or stop

YM

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