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

Last change on this file since 2616 was 1524, checked in by emillour, 9 years ago

All GCMS:
More updates to enforce dynamics/physics separation:

get rid of references to "temps_mod" from physics packages;
make a "time_phylmdz_mod.F90" module to store that
information and fill it via "iniphysiq".

EM

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