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

Last change on this file since 3026 was 3004, checked in by emillour, 18 months ago

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