module mo_schu !!$ use ppgrid, only : pverp !!$ use mo_grid, only : plevp use shr_kind_mod, only: r8 => shr_kind_r8 implicit none integer, parameter :: ngast = 17 integer, parameter :: tdim = 501 real(r8), parameter :: t_del = 5._r8/(tdim-1) real(r8), parameter :: t_fac = (tdim-1)/5._r8 integer :: ii, jj real(r8) :: d_table(0:tdim,ngast-1) real(r8) :: x_table(0:tdim,ngast-1) real(r8) :: o2_table(tdim) real(r8), dimension(12,ngast-1) :: a, b !----------------------------------------------------------------------------- ! a(16,12) coefficients for rj(m) (table 1 in kockarts 1994) ! b(16,12) rj(o2)(table 2 in kockarts 1994) ! rjm attenuation coefficients rj(m) ! rjo2 rj(o2) !----------------------------------------------------------------------------- data ((a(jj,ii),jj=1,12),ii=1,ngast-1) / & !a 57000-56500.5 cm-1 1.13402e-01_r8,1.00088e-20_r8,3.48747e-01_r8,2.76282e-20_r8,3.47322e-01_r8,1.01267e-19_r8, & 1.67351e-01_r8,5.63588e-19_r8,2.31433e-02_r8,1.68267e-18_r8,0.00000e+00_r8,0.00000e+00_r8, & !a 56500-56000.5 cm-1 2.55268e-03_r8,1.64489e-21_r8,1.85483e-01_r8,2.03591e-21_r8,2.60603e-01_r8,4.62276e-21_r8, & 2.50337e-01_r8,1.45106e-20_r8,1.92340e-01_r8,7.57381e-20_r8,1.06363e-01_r8,7.89634e-19_r8, & !a 56000-55500.5 cm-1 4.21594e-03_r8,8.46639e-22_r8,8.91886e-02_r8,1.12935e-21_r8,2.21334e-01_r8,1.67868e-21_r8, & 2.84446e-01_r8,3.94782e-21_r8,2.33442e-01_r8,1.91554e-20_r8,1.63433e-01_r8,2.25346e-19_r8, & !a 55500-55000.5 cm-1 3.93529e-03_r8,6.79660e-22_r8,4.46906e-02_r8,9.00358e-22_r8,1.33060e-01_r8,1.55952e-21_r8, & 3.25506e-01_r8,3.43763e-21_r8,2.79405e-01_r8,1.62086e-20_r8,2.10316e-01_r8,1.53883e-19_r8, & !a 55000-54500.5 cm-1 2.60939e-03_r8,2.33791e-22_r8,2.08101e-02_r8,3.21734e-22_r8,1.67186e-01_r8,5.77191e-22_r8, & 2.80694e-01_r8,1.33362e-21_r8,3.26867e-01_r8,6.10533e-21_r8,1.96539e-01_r8,7.83142e-20_r8, & !a 54500-54000.5 cm-1 9.33711e-03_r8,1.32897e-22_r8,3.63980e-02_r8,1.78786e-22_r8,1.46182e-01_r8,3.38285e-22_r8, & 3.81762e-01_r8,8.93773e-22_r8,2.58549e-01_r8,4.28115e-21_r8,1.64773e-01_r8,4.67537e-20_r8, & !a 54000-53500.5 cm-1 9.51799e-03_r8,1.00252e-22_r8,3.26320e-02_r8,1.33766e-22_r8,1.45962e-01_r8,2.64831e-22_r8, & 4.49823e-01_r8,6.42879e-22_r8,2.14207e-01_r8,3.19594e-21_r8,1.45616e-01_r8,2.77182e-20_r8, & !a 53500-53000.5 cm-1 7.87331e-03_r8,3.38291e-23_r8,6.91451e-02_r8,4.77708e-23_r8,1.29786e-01_r8,8.30805e-23_r8, & 3.05103e-01_r8,2.36167e-22_r8,3.35007e-01_r8,8.59109e-22_r8,1.49766e-01_r8,9.63516e-21_r8, & !a 53000-52500.5 cm-1 6.92175e-02_r8,1.56323e-23_r8,1.44403e-01_r8,3.03795e-23_r8,2.94489e-01_r8,1.13219e-22_r8, & 3.34773e-01_r8,3.48121e-22_r8,9.73632e-02_r8,2.10693e-21_r8,5.94308e-02_r8,1.26195e-20_r8, & !a 52500-52000.5 cm-1 1.47873e-01_r8,8.62033e-24_r8,3.15881e-01_r8,3.51859e-23_r8,4.08077e-01_r8,1.90524e-22_r8, & 8.08029e-02_r8,9.93062e-22_r8,3.90399e-02_r8,6.38738e-21_r8,8.13330e-03_r8,9.93644e-22_r8, & !a 52000-51500.5 cm-1 1.50269e-01_r8,1.02621e-23_r8,2.39823e-01_r8,3.48120e-23_r8,3.56408e-01_r8,1.69494e-22_r8, & 1.61277e-01_r8,6.59294e-22_r8,8.89713e-02_r8,2.94571e-21_r8,3.25063e-03_r8,1.25548e-20_r8, & !a 51500-51000.5 cm-1 2.55746e-01_r8,8.49877e-24_r8,2.94733e-01_r8,2.06878e-23_r8,2.86382e-01_r8,9.30992e-23_r8, & 1.21011e-01_r8,3.66239e-22_r8,4.21105e-02_r8,1.75700e-21_r8,0.00000e+00_r8,0.00000e+00_r8, & !a 51000-50500.5 cm-1 5.40111e-01_r8,7.36085e-24_r8,2.93263e-01_r8,2.46742e-23_r8,1.63417e-01_r8,1.37832e-22_r8, & 3.23781e-03_r8,2.15052e-21_r8,0.00000e+00_r8,0.00000e+00_r8,0.00000e+00_r8,0.00000e+00_r8, & !a 50500-50000.5 cm-1 8.18514e-01_r8,7.17937e-24_r8,1.82262e-01_r8,4.17496e-23_r8,0.00000e+00_r8,0.00000e+00_r8, & 0.00000e+00_r8,0.00000e+00_r8,0.00000e+00_r8,0.00000e+00_r8,0.00000e+00_r8,0.00000e+00_r8, & !a 50000-49500.5 cm-1 8.73680e-01_r8,7.13444e-24_r8,1.25583e-01_r8,2.77819e-23_r8,0.00000e+00_r8,0.00000e+00_r8, & 0.00000e+00_r8,0.00000e+00_r8,0.00000e+00_r8,0.00000e+00_r8,0.00000e+00_r8,0.00000e+00_r8, & !a 49500-49000.5 cm-1 3.32476e-04_r8,7.00362e-24_r8,9.89000e-01_r8,6.99600e-24_r8,0.00000e+00_r8,0.00000e+00_r8, & 0.00000e+00_r8,0.00000e+00_r8,0.00000e+00_r8,0.00000e+00_r8,0.00000e+00_r8,0.00000e+00_r8 / data ((b(jj,ii),jj=1,12),ii=1,ngast-1) / & ! 57000-56500.5 cm-1 1.07382e-21_r8,9.95029e-21_r8,7.19430e-21_r8,2.48960e-20_r8,2.53735e-20_r8,7.54467e-20_r8, & 4.48987e-20_r8,2.79981e-19_r8,9.72535e-20_r8,9.29745e-19_r8,2.30892e-20_r8,4.08009e-17_r8, & ! 56500-56000.5 cm-1 3.16903e-22_r8,1.98251e-21_r8,5.87326e-22_r8,3.44057e-21_r8,2.53094e-21_r8,8.81484e-21_r8, & 8.82299e-21_r8,4.17179e-20_r8,2.64703e-20_r8,2.43792e-19_r8,8.73831e-20_r8,1.46371e-18_r8, & ! 56000-55500.5 cm-1 1.64421e-23_r8,9.26011e-22_r8,2.73137e-22_r8,1.33640e-21_r8,9.79188e-22_r8,2.99706e-21_r8, & 3.37768e-21_r8,1.39438e-20_r8,1.47898e-20_r8,1.04322e-19_r8,4.08014e-20_r8,6.31023e-19_r8, & ! 55500-55000.5 cm-1 8.68729e-24_r8,7.31056e-22_r8,8.78313e-23_r8,1.07173e-21_r8,8.28170e-22_r8,2.54986e-21_r8, & 2.57643e-21_r8,9.42698e-21_r8,9.92377e-21_r8,5.21402e-20_r8,3.34301e-20_r8,2.91785e-19_r8, & ! 55000-54500.5 cm-1 1.20679e-24_r8,2.44092e-22_r8,2.64326e-23_r8,4.03998e-22_r8,2.53514e-22_r8,8.53166e-22_r8, & 1.29834e-21_r8,3.74482e-21_r8,5.12103e-21_r8,2.65798e-20_r8,2.10948e-20_r8,2.35315e-19_r8, & ! 54500-54000.5 cm-1 2.79656e-24_r8,1.40820e-22_r8,3.60824e-23_r8,2.69510e-22_r8,4.02850e-22_r8,8.83735e-22_r8, & 1.77198e-21_r8,6.60221e-21_r8,9.60992e-21_r8,8.13558e-20_r8,4.95591e-21_r8,1.22858e-17_r8, & ! 54000-53500.5 cm-1 2.36959e-24_r8,1.07535e-22_r8,2.83333e-23_r8,2.16789e-22_r8,3.35242e-22_r8,6.42753e-22_r8, & 1.26395e-21_r8,5.43183e-21_r8,4.88083e-21_r8,5.42670e-20_r8,3.27481e-21_r8,1.58264e-17_r8, & ! 53500-53000.5 cm-1 8.65018e-25_r8,3.70310e-23_r8,1.04351e-23_r8,6.43574e-23_r8,1.17431e-22_r8,2.70904e-22_r8, & 4.88705e-22_r8,1.65505e-21_r8,2.19776e-21_r8,2.71172e-20_r8,2.65257e-21_r8,2.13945e-17_r8, & ! 53000-52500.5 cm-1 9.63263e-25_r8,1.54249e-23_r8,4.78065e-24_r8,2.97642e-23_r8,6.40637e-23_r8,1.46464e-22_r8, & 1.82634e-22_r8,7.12786e-22_r8,1.64805e-21_r8,2.37376e-17_r8,9.33059e-22_r8,1.13741e-20_r8, & ! 52500-52000.5 cm-1 1.08414e-24_r8,8.37560e-24_r8,9.15550e-24_r8,2.99295e-23_r8,9.38405e-23_r8,1.95845e-22_r8, & 2.84356e-22_r8,3.39699e-21_r8,1.94524e-22_r8,2.72227e-19_r8,1.18924e-21_r8,3.20246e-17_r8, & ! 52000-51500.5 cm-1 1.52817e-24_r8,1.01885e-23_r8,1.22946e-23_r8,4.16517e-23_r8,9.01287e-23_r8,2.34869e-22_r8, & 1.93510e-22_r8,1.44956e-21_r8,1.81051e-22_r8,5.17773e-21_r8,9.82059e-22_r8,6.22768e-17_r8, & ! 51500-51000.5 cm-1 2.12813e-24_r8,8.48035e-24_r8,5.23338e-24_r8,1.93052e-23_r8,1.99464e-23_r8,7.48997e-23_r8, & 4.96642e-22_r8,6.15691e-17_r8,4.47504e-23_r8,2.76004e-22_r8,8.26788e-23_r8,1.65278e-21_r8, & ! 51000-50500.5 cm-1 3.81336e-24_r8,7.32307e-24_r8,5.60549e-24_r8,2.04651e-23_r8,3.36883e-22_r8,6.15708e-17_r8, & 2.09877e-23_r8,1.07474e-22_r8,9.13562e-24_r8,8.41252e-22_r8,0.00000e+00_r8,0.00000e+00_r8, & ! 50500-50000.5 cm-1 5.75373e-24_r8,7.15986e-24_r8,5.90031e-24_r8,3.05375e-23_r8,2.97196e-22_r8,8.92000e-17_r8, & 8.55920e-24_r8,1.66709e-17_r8,0.00000e+00_r8,0.00000e+00_r8,0.00000e+00_r8,0.00000e+00_r8, & ! 50000-49500.5 cm-1 6.21281e-24_r8,7.13108e-24_r8,3.30780e-24_r8,2.61196e-23_r8,1.30783e-22_r8,9.42550e-17_r8, & 2.69241e-24_r8,1.46500e-17_r8,0.00000e+00_r8,0.00000e+00_r8,0.00000e+00_r8,0.00000e+00_r8, & ! 49500-49000.5 cm-1 6.81118e-24_r8,6.98767e-24_r8,7.55667e-25_r8,2.75124e-23_r8,1.94044e-22_r8,1.45019e-16_r8, & 1.92236e-24_r8,3.73223e-17_r8,0.00000e+00_r8,0.00000e+00_r8,0.00000e+00_r8,0.00000e+00_r8 / contains subroutine schu_inti !----------------------------------------------------------------------------- ! ... initialize the tables !----------------------------------------------------------------------------- implicit none !----------------------------------------------------------------------------- ! ... local variables !----------------------------------------------------------------------------- integer :: i, iw, k, j1, jp1, j real(r8) :: col real(r8), dimension(6) :: a0, a1, b0, b1 do iw = 1,ngast-1 x_table(0,iw) = sum( a(1:11:2,iw) ) d_table(0,iw) = sum( b(1:11:2,iw) ) do k = 1,tdim col = 22._r8 + t_del*real(k-1) o2_table(k) = col col = 10._r8**col a1(:) = a(2:12:2,iw)*col b1(:) = b(2:12:2,iw)*col where( a1(:) < 500._r8 ) a0(:) = exp( -a1(:) ) elsewhere a0(:) = 0._r8 endwhere where( b1(:) < 500._r8 ) b0(:) = exp( -b1(:) ) elsewhere b0(:) = 0._r8 endwhere x_table(k,iw) = dot_product( a(1:11:2,iw),a0(:) ) d_table(k,iw) = dot_product( b(1:11:2,iw),b0(:) ) end do end do end subroutine schu_inti subroutine schu( o2col, secchi, dto2, xscho2 ) !----------------------------------------------------------------------------- ! purpose: ! calculate the equivalent absorption cross section of o2 in the sr bands. ! the algorithm is based on: g.kockarts, penetration of solar radiation ! in the schumann-runge bands of molecular oxygen: a robust approximation, ! annales geophysicae, v12, n12, pp. 1207ff, dec 1994. calculation is ! done on the wavelength grid used by kockarts (1994). final values do ! include effects from the herzberg continuum. !----------------------------------------------------------------------------- ! parameters: ! nz - integer, number of specified altitude levels in the working (i) ! grid ! o2col - real, slant overhead o2 column (molec/cc) at each specified (i) ! altitude ! dto2 - real, optical depth due to o2 absorption at each specified (o) ! vertical layer at each specified wavelength ! xscho2 - real, molecular absorption cross section in sr bands at (o) ! each specified wavelength. includes herzberg continuum !----------------------------------------------------------------------------- !!$ use mo_grid, only : plev, plevp use ppgrid, only : pver,pverp implicit none !----------------------------------------------------------------------------- ! ... dummy arguments !----------------------------------------------------------------------------- real(r8), intent(inout) :: dto2(pver,ngast-1) real(r8), intent(inout) :: xscho2(pverp,ngast-1) real(r8), intent(in) :: o2col(pverp) real(r8), intent(in) :: secchi(pverp) !----------------------------------------------------------------------------- ! ... local variables !----------------------------------------------------------------------------- integer :: i, iw, j, j1, jp1, k, ki, kp1 integer :: index(pverp) integer :: minind(1), maxind(1) real(r8) :: a0, a1, b0, b1 real(r8), dimension(6) :: ac, bc real(r8), dimension(pverp,6) :: aa, bb real(r8), dimension(pverp) :: rjm, rjo2 real(r8), dimension(pverp) :: rjmi, rjo2i, lo2col, dels !----------------------------------------------------------------------------- ! ... initialize r(m) !----------------------------------------------------------------------------- rjm(1:pverp) = 0._r8 rjo2(1:pverp) = 0._r8 !----------------------------------------------------------------------------- ! ... initialize the table interpolation !----------------------------------------------------------------------------- where( o2col(:) /= 0 ) lo2col(:) = log10( o2col(:) ) endwhere do ki = 1,pverp if( o2col(ki) /= 0._r8 ) then if( lo2col(ki) <= o2_table(1) ) then dels(ki) = 0._r8 index(ki) = 1 else if( lo2col(ki) >= o2_table(tdim) ) then dels(ki) = 1._r8 index(ki) = tdim-1 else do k = 2,tdim if( lo2col(ki) <= o2_table(k) ) then dels(ki) = t_fac*(lo2col(ki) - o2_table(k-1)) index(ki) = k-1 exit end if end do end if else index(ki) = 0 dels(ki) = 0._r8 end if end do !----------------------------------------------------------------------------- ! ... calculate sum of exponentials (eqs 7 and 8 of kockarts 1994) !----------------------------------------------------------------------------- do iw = 1,ngast-1 do k = 1,pverp ki = index(k) rjm(k) = x_table(ki,iw) + dels(k)*(x_table(ki+1,iw) - x_table(ki,iw)) rjo2(k) = d_table(ki,iw) + dels(k)*(d_table(ki+1,iw) - d_table(ki,iw)) end do do k = 1,pver if( rjm(k) > 1.e-100_r8 ) then kp1 = k + 1 if( rjm(kp1) > 0._r8 ) then dto2(k,iw) = log( rjm(kp1) ) / secchi(kp1) - log( rjm(k) ) * secchi(k) else dto2(k,iw) = 1000._r8 end if else dto2(k,iw) = 1000._r8 end if end do do k = 1,pverp if( rjm(k) > 1.e-100_r8 ) then if( rjo2(k) > 1.e-100_r8 ) then xscho2(k,iw) = rjo2(k)/rjm(k) else xscho2(k,iw) = 0._r8 end if else xscho2(k,iw) = 0._r8 end if end do end do end subroutine schu end module mo_schu