| 1 | module lwb_mod |
|---|
| 2 | |
|---|
| 3 | implicit none |
|---|
| 4 | |
|---|
| 5 | contains |
|---|
| 6 | |
|---|
| 7 | subroutine lwb (kdlon,kflev,tlev,tlay,dt0 |
|---|
| 8 | . ,bsurf,btop,blay,blev,dblay,dbsublay) |
|---|
| 9 | |
|---|
| 10 | c---------------------------------------------------------------------- |
|---|
| 11 | c LWB computes the planck function and gradient |
|---|
| 12 | c from a polynomial development of planck function |
|---|
| 13 | c---------------------------------------------------------------------- |
|---|
| 14 | |
|---|
| 15 | use dimradmars_mod, only: ndlon, ndlo2, nir |
|---|
| 16 | use yomlw_h, only: nlaylte, xi , tstand, xp |
|---|
| 17 | implicit none |
|---|
| 18 | |
|---|
| 19 | c---------------------------------------------------------------------- |
|---|
| 20 | c 0.1 arguments |
|---|
| 21 | c --------- |
|---|
| 22 | c inputs: |
|---|
| 23 | c ------- |
|---|
| 24 | integer kdlon ! part of ngrid |
|---|
| 25 | integer kflev ! part of nalyer |
|---|
| 26 | |
|---|
| 27 | real dt0 (ndlo2) ! surface temperature discontinuity |
|---|
| 28 | real tlay (ndlo2,kflev) ! layer temperature |
|---|
| 29 | real tlev (ndlo2,kflev+1) ! level temperature |
|---|
| 30 | |
|---|
| 31 | c outputs: |
|---|
| 32 | c -------- |
|---|
| 33 | real bsurf (ndlo2,nir) ! surface spectral planck function |
|---|
| 34 | real btop (ndlo2,nir) ! top spectral planck function |
|---|
| 35 | real blev (ndlo2,nir,kflev+1) ! level spectral planck function |
|---|
| 36 | real blay (ndlo2,nir,kflev) ! layer spectral planck function |
|---|
| 37 | real dblay (ndlo2,nir,kflev) ! layer gradient spectral planck function |
|---|
| 38 | real dbsublay (ndlo2,nir,2*kflev) ! layer gradient spectral planck |
|---|
| 39 | ! function in sub layers |
|---|
| 40 | |
|---|
| 41 | c---------------------------------------------------------------------- |
|---|
| 42 | c 0.2 local arrays |
|---|
| 43 | c ------------ |
|---|
| 44 | |
|---|
| 45 | integer jk, jl, jnu, jk1, jk2 |
|---|
| 46 | |
|---|
| 47 | real ztlay (ndlon) |
|---|
| 48 | real ztlev (ndlon) |
|---|
| 49 | |
|---|
| 50 | c---------------------------------------------------------------------- |
|---|
| 51 | do jnu=1,nir |
|---|
| 52 | c---------------------------------------------------------------------- |
|---|
| 53 | c 1.1 levels and layers from surface to nlaylte |
|---|
| 54 | c --------------------------------------- |
|---|
| 55 | |
|---|
| 56 | do jk = 1 , nlaylte |
|---|
| 57 | do jl = 1 , kdlon |
|---|
| 58 | |
|---|
| 59 | c level planck function |
|---|
| 60 | c --------------------- |
|---|
| 61 | ztlev(jl)=(tlev(jl,jk)-tstand)/tstand ! tstand = 200k |
|---|
| 62 | |
|---|
| 63 | blev(jl,jnu,jk) = xp(1,jnu) |
|---|
| 64 | . +ztlev(jl)*(xp(2,jnu) |
|---|
| 65 | . +ztlev(jl)*(xp(3,jnu) |
|---|
| 66 | . +ztlev(jl)*(xp(4,jnu) |
|---|
| 67 | . +ztlev(jl)*(xp(5,jnu) |
|---|
| 68 | . +ztlev(jl)*(xp(6,jnu) ))))) |
|---|
| 69 | |
|---|
| 70 | c layer planck function |
|---|
| 71 | c --------------------- |
|---|
| 72 | ztlay(jl)=(tlay(jl,jk)-tstand)/tstand |
|---|
| 73 | |
|---|
| 74 | blay(jl,jnu,jk) = xp(1,jnu) |
|---|
| 75 | . +ztlay(jl)*(xp(2,jnu) |
|---|
| 76 | . +ztlay(jl)*(xp(3,jnu) |
|---|
| 77 | . +ztlay(jl)*(xp(4,jnu) |
|---|
| 78 | . +ztlay(jl)*(xp(5,jnu) |
|---|
| 79 | . +ztlay(jl)*(xp(6,jnu) ))))) |
|---|
| 80 | |
|---|
| 81 | c planck function gradient |
|---|
| 82 | c ------------------------ |
|---|
| 83 | dblay(jl,jnu,jk) = xp(2,jnu) |
|---|
| 84 | . +ztlay(jl)*(2*xp(3,jnu) |
|---|
| 85 | . +ztlay(jl)*(3*xp(4,jnu) |
|---|
| 86 | . +ztlay(jl)*(4*xp(5,jnu) |
|---|
| 87 | . +ztlay(jl)*(5*xp(6,jnu) )))) |
|---|
| 88 | dblay(jl,jnu,jk) = dblay(jl,jnu,jk)/tstand |
|---|
| 89 | |
|---|
| 90 | enddo |
|---|
| 91 | enddo |
|---|
| 92 | |
|---|
| 93 | c---------------------------------------------------------------------- |
|---|
| 94 | c 1.2 top of the atmosphere and surface |
|---|
| 95 | c -------------------------------- |
|---|
| 96 | |
|---|
| 97 | do jl = 1 , kdlon |
|---|
| 98 | c top of the atmosphere |
|---|
| 99 | c --------------------- |
|---|
| 100 | ztlev(jl) = (tlev(jl,nlaylte+1)-tstand)/tstand |
|---|
| 101 | |
|---|
| 102 | blev(jl,jnu,nlaylte+1) = xp(1,jnu) |
|---|
| 103 | . +ztlev(jl)*(xp(2,jnu) |
|---|
| 104 | . +ztlev(jl)*(xp(3,jnu) |
|---|
| 105 | . +ztlev(jl)*(xp(4,jnu) |
|---|
| 106 | . +ztlev(jl)*(xp(5,jnu) |
|---|
| 107 | . +ztlev(jl)*(xp(6,jnu) ))))) |
|---|
| 108 | btop(jl,jnu) = blev(jl,jnu,nlaylte+1) |
|---|
| 109 | |
|---|
| 110 | c surface |
|---|
| 111 | c ------- |
|---|
| 112 | ztlay(jl) = (tlev(jl,1)+dt0(jl)-tstand)/tstand |
|---|
| 113 | |
|---|
| 114 | bsurf(jl,jnu) = xp(1,jnu) |
|---|
| 115 | . +ztlay(jl)*(xp(2,jnu) |
|---|
| 116 | . +ztlay(jl)*(xp(3,jnu) |
|---|
| 117 | . +ztlay(jl)*(xp(4,jnu) |
|---|
| 118 | . +ztlay(jl)*(xp(5,jnu) |
|---|
| 119 | . +ztlay(jl)*(xp(6,jnu) ))))) |
|---|
| 120 | |
|---|
| 121 | enddo |
|---|
| 122 | |
|---|
| 123 | c---------------------------------------------------------------------- |
|---|
| 124 | c 1.3 Gradients in sub-layers |
|---|
| 125 | c ----------------------- |
|---|
| 126 | |
|---|
| 127 | do jk=1,nlaylte |
|---|
| 128 | jk2 = 2 * jk |
|---|
| 129 | jk1 = jk2 - 1 |
|---|
| 130 | do jl=1,kdlon |
|---|
| 131 | dbsublay(jl,jnu,jk1)=blay(jl,jnu,jk)-blev(jl,jnu,jk) |
|---|
| 132 | dbsublay(jl,jnu,jk2)=blev(jl,jnu,jk+1)-blay(jl,jnu,jk) |
|---|
| 133 | enddo |
|---|
| 134 | enddo |
|---|
| 135 | |
|---|
| 136 | c---------------------------------------------------------------------- |
|---|
| 137 | enddo ! (do jnu=1,nir) |
|---|
| 138 | c---------------------------------------------------------------------- |
|---|
| 139 | |
|---|
| 140 | end subroutine lwb |
|---|
| 141 | |
|---|
| 142 | end module lwb_mod |
|---|