source: trunk/mars/libf/phymars/lwb.F @ 38

Last change on this file since 38 was 38, checked in by emillour, 14 years ago

Ajout du modè Martien (mon LMDZ.MARS.BETA, du 28/01/2011) dans le rértoire mars, pour pouvoir suivre plus facilement les modifs.
EM

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