source: trunk/LMDZ.MARS/libf/phymars/lwi.F @ 4022

Last change on this file since 4022 was 3726, checked in by emillour, 11 months ago

Mars PCM:
Turn "callkeys.h" into module "callkeys_mod.F90"
EM

File size: 5.5 KB
Line 
1      module lwi_mod
2     
3      implicit none
4     
5      contains
6
7      subroutine lwi (ig0,kdlon,kflev
8     .                ,psi,zdblay,pdp
9     .                ,newpcolc )
10
11
12      use dimradmars_mod, only: ndlo2, ndlon, nflev, nir
13      use yomlw_h, only: gcp, nlaylte, xi
14      USE comcstfi_h, ONLY: g, cpp
15      USE time_phylmdz_mod, ONLY: dtphys
16      use callkeys_mod, only: semi, iradia
17      implicit none
18
19CCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCC
20C
21C                             -   lwi    -   
22C
23C     PURPOSE:       Shema semi - implicite
24C
25C
26CCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCC
27
28
29c************************************************************************
30c
31c        0.    Declarations
32c              ------------
33c
34c-------------------------------------------------------------------------
35c        0.1   Arguments
36c              ---------
37c
38 
39      integer,intent(in) :: ig0
40      integer,intent(in) :: kdlon
41      integer,intent(in) :: kflev
42
43      real,intent(in) :: psi(ndlo2,kflev)
44      real,intent(in) :: zdblay(ndlo2,nir,kflev)
45      real,intent(in) :: pdp(ndlo2,kflev)
46
47      real,intent(out) :: newpcolc(ndlo2,kflev)
48
49c-------------------------------------------------------------------------
50c        0.2   local arrays
51c              ------------
52c
53      real    di(ndlon,nflev)
54     .      , hi(ndlon,nflev)
55     .      , bi(ndlon,nflev)
56
57      real    ci(ndlon,nflev)
58     .      , ai(ndlon,nflev)
59      real    deltat
60
61      real   semit, denom
62
63      integer i, jl
64
65c************************************************************************
66c
67c        1.    Initialisations
68c              ---------------
69c
70c-----------------------------------------------------------------------
71 
72        deltat = dtphys * iradia
73c       print*,'SEMI = ',semi, '(expl:0  semi-implicite:0.5  impl:1)'
74        semit = semi * deltat
75c       semi = 0.
76
77c       print*,'dtphys,iradia,deltat,semit:',dtphys,iradia,deltat,semit
78c       print*,'g,cpp',g,cpp
79
80
81c************************************************************************
82c
83c        2.   
84c              ---------------
85c
86c-------------------------------------------------------------------------
87c        2.1   Calcul des di
88c              -------------
89c
90
91
92      do i = 1 , nlaylte-1
93        do jl = 1 , kdlon
94c     -------------------
95      di(jl,i) =  1 + semit * (g / pdp(jl,i) / cpp) * (
96     .    ( xi(ig0+jl,1,i,nlaylte+1)
97     .    + xi(ig0+jl,1,i,i+1)
98     .    + xi(ig0+jl,1,i,i-1) )
99     .    *    zdblay(jl,1,i)
100     .  + ( xi(ig0+jl,2,i,nlaylte+1)
101     .    + xi(ig0+jl,2,i,i+1)
102     .    + xi(ig0+jl,2,i,i-1) )
103     .    *    zdblay(jl,2,i)
104     .     )
105c     -------------------
106        enddo
107      enddo
108
109c couche nlaylte
110c ------------
111c      , on enleve i,i+1 sinon on a 2 fois le cooling2space
112
113      do jl = 1 , kdlon
114c     -------------------
115      di(jl,nlaylte) =  1 + semit * (g / pdp(jl,nlaylte) / cpp) * (
116     .    ( xi(ig0+jl,1,nlaylte,nlaylte+1)
117     .    + xi(ig0+jl,1,nlaylte,nlaylte-1) )
118     .    *    zdblay(jl,1,nlaylte)
119     .  + ( xi(ig0+jl,2,nlaylte,nlaylte+1)
120     .    + xi(ig0+jl,2,nlaylte,nlaylte-1) )
121     .    *    zdblay(jl,2,nlaylte)
122     .     )
123c     -------------------
124      enddo
125
126c-------------------------------------------------------------------------
127c        2.2   Calcul des hi
128c              -------------
129c
130
131      do i = 1 , nlaylte-1
132        do jl = 1 , kdlon
133c     -------------------
134      hi(jl,i) =    - semit * (g / pdp(jl,i) / cpp) *
135     .            (    xi(ig0+jl,1,i,i+1) * zdblay(jl,1,i+1)   
136     .               + xi(ig0+jl,2,i,i+1) * zdblay(jl,2,i+1)   )
137c     -------------------
138        enddo
139      enddo
140
141c-------------------------------------------------------------------------
142c        2.3   Calcul des bi
143c              -------------
144c
145
146
147      do i = 2 , nlaylte
148        do jl = 1 , kdlon
149c     -------------------
150      bi(jl,i) =   - semit * (g / pdp(jl,i) / cpp) *
151     .           (     xi(ig0+jl,1,i,i-1) * zdblay(jl,1,i-1)   
152     .               + xi(ig0+jl,2,i,i-1) * zdblay(jl,2,i-1)   )
153c     -------------------
154        enddo
155      enddo
156
157
158c couche 1
159c --------
160c  tant qu'on a pas un calcul propre de zdblay(0) qui tienne compte de
161c    la discontinuite de temperature au sol , on met  b1 = 0
162
163
164      do jl = 1 , kdlon
165        bi(jl,1) = 0
166      enddo
167
168c-------------------------------------------------------------------------
169c        2.4   
170c              -------------
171c
172
173c couche nlaylte
174c ------------
175
176      do jl = 1 , kdlon
177c     -------------------
178      ci(jl,nlaylte) = (gcp * psi(jl,nlaylte) / pdp(jl,nlaylte))
179     .                   / di(jl,nlaylte)
180
181      ai(jl,nlaylte) = - bi(jl,nlaylte) / di(jl,nlaylte)
182c     -------------------
183      enddo
184
185
186
187      do i = nlaylte-1 , 1 , -1
188        do jl = 1 , kdlon
189c     -------------------
190      denom = di(jl,i) + hi(jl,i) * ai(jl,i+1)
191
192      ci(jl,i) = (  gcp * psi(jl,i) / pdp(jl,i)
193     .             - hi(jl,i) * ci(jl,i+1)  )  / denom
194 
195      ai(jl,i) = -bi(jl,i) / denom
196c     -------------------
197        enddo
198      enddo
199
200
201c-------------------------------------------------------------------------
202c        2.5   
203c              -------------
204c
205
206c couche 1
207c -------
208      do jl = 1 , kdlon
209        newpcolc(jl,1) = ci(jl,1)
210      enddo
211
212
213      do i = 2 , nlaylte
214        do jl = 1 , kdlon
215           newpcolc(jl,i) = ci(jl,i) + ai(jl,i) * newpcolc(jl,i-1)
216        enddo
217      enddo
218
219
220
221c-------------------------------------------------------------------------
222
223      end subroutine lwi
224     
225      end module lwi_mod
Note: See TracBrowser for help on using the repository browser.