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

Last change on this file since 3026 was 1266, checked in by aslmd, 11 years ago

LMDZ.MARS
IMPORTANT CHANGE

  • Remove all reference/use of nlayermx and dimphys.h
  • Made use of automatic arrays whenever arrays are needed with dimension nlayer
  • Remove lots of obsolete reference to dimensions.h
  • Converted iono.h and param_v4.h into corresponding modules

(with embedded subroutine to allocate arrays)
(no arrays allocated if thermosphere not used)

  • Deleted param.h and put contents into module param_v4_h
  • Adapted testphys1d, newstart, etc...
  • Made DATA arrays in param_read to be initialized by subroutine

fill_data_thermos in module param_v4_h

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