source: trunk/LMDZ.PLUTO/libf/phypluto/setspi.F90 @ 3504

Last change on this file since 3504 was 3504, checked in by afalco, 2 weeks ago

Pluto: print only on master process.
AF

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