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