source: trunk/LMDZ.MARS/libf/phymars/lwb.F @ 1085

Last change on this file since 1085 was 1085, checked in by emillour, 11 years ago

Removed wrongly referenced module variables in lwb and physdem.
Added option "-full" to makegcm_* to force a full recompilation from scratch (ie: removing libo/... and makefile before launching makefile building and compilation).
EM

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