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

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

LMDZ.MARS. related to r1266. forgot to remove a few now-obsolete dimensions.h includes in Mars physics.

File size: 13.5 KB
Line 
1      SUBROUTINE flusv(KDLON,nsf,n,omega,g,tau,emis,bh,bsol,fah,fdh)
2      use dimradmars_mod, only: ndlo2, ndlon, nflev
3      IMPLICIT NONE
4c.......................................................................
5c
6c  calcul des flux ascendant et descendant aux interfaces entre n couches
7c  * dans l'infrarouge
8c  * B est une fonction lineaire de $\tau$ a l'interieur de chaque couche
9c  * le B du sol peut etre different de celui qui correspond au profil
10c    de la n-ieme couche
11c  * l'hypothese est une hypothese a deux flux isotropes sur chaque
12c    hemisphere ("hemispheric constant") + "source function technique"
13c    (voir Toon et al. 1988)
14c  * le flux descendant en haut de l'atmosphere est nul
15c  * les couches sont numerotees du haut de l'atmosphere vers le sol
16c
17c  in :   * KDLON      ---> dimension de vectorisation
18c         * nsf        ---> nsf=0 ==> "hemispheric constant"
19c                           nsf>0 ==> "hemispheric constant" + "source function"
20c         * n          ---> nombre de couches
21c         * omega(i)   ---> single scattering albedo pour la i-eme couche
22c         * g(i)       ---> asymmetry parameter pour la i-eme couche
23c         * tau(i)     ---> epaisseur optique de la i-eme couche
24c         * emis       ---> emissivite du sol
25c         * bh(i)      ---> luminance du corps noir en haut de la i-eme
26c                           couche, bh(n+1) pour la valeur au sol qui
27c                           correspond au profil de la n-ieme couche
28c         * bsol       ---> luminance du corps noir au sol
29c
30c  out :  * fah(i)     ---> flux ascendant en haut de la i-eme couche,
31c                           fah(n+1) pour le sol
32c         * fdh(i)     ---> flux descendant en haut de la i-eme couche,
33c                           fdh(n+1) pour le sol
34c
35c.......................................................................
36c  declaration des arguments
37c
38      INTEGER KDLON,nsf,n
39      REAL omega(NDLO2,n),g(NDLO2,n),tau(NDLO2,n),emis(NDLO2)
40     &,bh(NDLO2,n+1),bsol(NDLO2),fah(NDLO2,n+1),fdh(NDLO2,n+1)
41c.......................................................................
42c  declaration des variables locales
43c
44      REAL pi
45      PARAMETER (pi=3.141592653589793E+0)
46      INTEGER iv,i,j
47      REAL beta,gama1,gama2,amu1,grgama,b0,b1
48      REAL a(NDLON,4*nflev),b(NDLON,4*nflev)
49     &    ,d(NDLON,4*nflev),e(NDLON,4*nflev)
50     &    ,y(NDLON,4*nflev)
51     &    ,alambda(NDLON,2*nflev)
52     &    ,e1(NDLON,2*nflev),e2(NDLON,2*nflev)
53     &    ,e3(NDLON,2*nflev),e4(NDLON,2*nflev)
54     &    ,cah(NDLON,2*nflev),cab(NDLON,2*nflev)
55     &    ,cdh(NDLON,2*nflev),cdb(NDLON,2*nflev)
56      REAL grg(NDLON,2*nflev),grh(NDLON,2*nflev)
57     &    ,grj(NDLON,2*nflev),grk(NDLON,2*nflev)
58     &    ,alpha1(NDLON,2*nflev),alpha2(NDLON,2*nflev)
59     &    ,sigma1(NDLON,2*nflev),sigma2(NDLON,2*nflev)
60      INTEGER nq
61      PARAMETER (nq=8)
62      REAL x(nq),w(nq),gri(NDLON,nq)
63      DATA x/1.9855071751231860E-2 , 0.1016667612931866E+0
64     &     , 0.2372337950418355E+0 , 0.4082826787521751E+0
65     &     , 0.5917173212478250E+0 , 0.7627662049581645E+0
66     &     , 0.8983332387068134E+0 , 0.9801449282487682E+0/
67      DATA w/5.0614268145185310E-2 , 0.1111905172266872E+0
68     &     , 0.1568533229389437E+0 , 0.1813418916891810E+0
69     &     , 0.1813418916891810E+0 , 0.1568533229389437E+0
70     &     , 0.1111905172266872E+0 , 5.0614268145185310E-2/
71c.......................................................................
72c
73c.......................................................................
74      do 10001 i=1,n
75                                                   do 99999 iv=1,KDLON
76      beta=(1.E+0-g(iv,i))/2.E+0
77      gama1=2.E+0*(1.E+0-omega(iv,i)*(1.E+0-beta))
78      gama2=2.E+0*omega(iv,i)*beta
79      amu1=5.E-1
80      alambda(iv,i)=sqrt(gama1**2-gama2**2)
81      grgama=(gama1-alambda(iv,i))/gama2
82c
83c ici une petite bidouille : si l'epaisseur optique d'une couche
84c est trop faible, $dB \over d\tau$ devient tres grand et le schema
85c s'ecroule. dans ces cas, on fait l'hypothese de couche isotherme.
86c
87      if (tau(iv,i).gt.1.E-3) then
88      b0=bh(iv,i)
89      b1=(bh(iv,i+1)-b0)/tau(iv,i)
90      else
91      b0=(bh(iv,i)+bh(iv,i+1))/2.E+0
92      b1=0.E+0
93      endif
94c
95      e1(iv,i)=1.E+0+grgama*exp(-alambda(iv,i)*tau(iv,i))
96      e2(iv,i)=1.E+0-grgama*exp(-alambda(iv,i)*tau(iv,i))
97      e3(iv,i)=grgama+exp(-alambda(iv,i)*tau(iv,i))
98      e4(iv,i)=grgama-exp(-alambda(iv,i)*tau(iv,i))
99      cah(iv,i)=2.E+0*pi*amu1*(b0+b1/(gama1+gama2))
100      cab(iv,i)=2.E+0*pi*amu1*(b0+b1*(tau(iv,i)+1.E+0/(gama1+gama2)))
101      cdh(iv,i)=2.E+0*pi*amu1*(b0-b1/(gama1+gama2))
102      cdb(iv,i)=2.E+0*pi*amu1*(b0+b1*(tau(iv,i)-1.E+0/(gama1+gama2))) 
103c
104      grg(iv,i)=(1.E+0/amu1-alambda(iv,i))
105      grh(iv,i)=grgama*(alambda(iv,i)+1.E+0/amu1)
106      grj(iv,i)=grh(iv,i)
107      grk(iv,i)=grg(iv,i)
108      alpha1(iv,i)=2.E+0*pi*(b0+b1*(1.E+0/(gama1+gama2)-amu1))
109      alpha2(iv,i)=2.E+0*pi*b1
110      sigma1(iv,i)=2.E+0*pi*(b0-b1*(1.E+0/(gama1+gama2)-amu1))
111      sigma2(iv,i)=alpha2(iv,i)
112c
11399999                                              continue
11410001 continue
115c.......................................................................
116                                                   do 99998 iv=1,KDLON
117      a(iv,1)=0.E+0
118      b(iv,1)=e1(iv,1)
119      d(iv,1)=-e2(iv,1)
120      e(iv,1)=-cdh(iv,1)
12199998                                              continue
122c
123      do 10002 i=1,n-1
124      j=2*i+1
125                                                   do 99997 iv=1,KDLON
126      a(iv,j)=e2(iv,i)*e3(iv,i)-e4(iv,i)*e1(iv,i)
127      b(iv,j)=e1(iv,i)*e1(iv,i+1)-e3(iv,i)*e3(iv,i+1)
128      d(iv,j)=e3(iv,i)*e4(iv,i+1)-e1(iv,i)*e2(iv,i+1)
129      e(iv,j)=e3(iv,i)*(cah(iv,i+1)-cab(iv,i))
130     &       +e1(iv,i)*(cdb(iv,i)-cdh(iv,i+1))
13199997                                              continue
13210002 continue
133c
134      do 10003 i=1,n-1
135      j=2*i
136                                                   do 99996 iv=1,KDLON
137      a(iv,j)=e2(iv,i+1)*e1(iv,i)-e3(iv,i)*e4(iv,i+1)
138      b(iv,j)=e2(iv,i)*e2(iv,i+1)-e4(iv,i)*e4(iv,i+1)
139      d(iv,j)=e1(iv,i+1)*e4(iv,i+1)-e2(iv,i+1)*e3(iv,i+1)
140      e(iv,j)=e2(iv,i+1)*(cah(iv,i+1)-cab(iv,i))
141     &       +e4(iv,i+1)*(cdb(iv,i)-cdh(iv,i+1))
14299996                                              continue
14310003 continue
144c     
145      j=2*n
146                                                   do 99995 iv=1,KDLON
147      a(iv,j)=e1(iv,n)-(1.E+0-emis(iv))*e3(iv,n)
148      b(iv,j)=e2(iv,n)-(1.E+0-emis(iv))*e4(iv,n)
149      d(iv,j)=0.E+0
150      e(iv,j)=emis(iv)*pi*bsol(iv)-cab(iv,n)+(1.E+0-emis(iv))*cdb(iv,n)
15199995                                              continue
152c.......................................................................
153      call sys3v(KDLON,2*n,a,b,d,e,y)
154c.......................................................................
155      do 10004 i=1,n
156                                                   do 99994 iv=1,KDLON
157      grg(iv,i)=grg(iv,i)*(y(iv,2*i-1)+y(iv,2*i))
158      grh(iv,i)=grh(iv,i)*(y(iv,2*i-1)-y(iv,2*i))
159      grj(iv,i)=grj(iv,i)*(y(iv,2*i-1)+y(iv,2*i))
160      grk(iv,i)=grk(iv,i)*(y(iv,2*i-1)-y(iv,2*i))
16199994                                              continue
16210004 continue
163c.......................................................................
164c les valeurs de flux "hemispheric constant"
165c
166      IF (nsf.eq.0) THEN
167      do 10005 i=1,n
168                                                   do 99993 iv=1,KDLON
169      fah(iv,i)=e3(iv,i)*y(iv,2*i-1)-e4(iv,i)*y(iv,2*i)+cah(iv,i)
170      fdh(iv,i)=e1(iv,i)*y(iv,2*i-1)-e2(iv,i)*y(iv,2*i)+cdh(iv,i)
17199993                                              continue
17210005 continue
173                                                   do 99992 iv=1,KDLON
174      fah(iv,n+1)=e1(iv,n)*y(iv,2*n-1)+e2(iv,n)*y(iv,2*n)+cab(iv,n)
175      fdh(iv,n+1)=e3(iv,n)*y(iv,2*n-1)+e4(iv,n)*y(iv,2*n)+cdb(iv,n)
17699992                                              continue
177      GOTO 10100
178      ENDIF
179c.......................................................................
180c passage a "source function"
181c
182c on applique une quadrature de dimension nq
183c x est le vecteur des \mu de la quadrature
184c w est le vecteur des poids correspondants
185c nq est fixe en parameter
186c x et w sont fixes en data
187c
188c.......................................................................
189c on part d'en haut et on descent selon les nq angles pour calculer
190c tous les flux descendants
191c
192      do 10006 j=1,nq
193                                                   do 99991 iv=1,KDLON
194      gri(iv,j)=0.E+0
19599991                                              continue
19610006 continue
197                                                   do 99990 iv=1,KDLON
198      fdh(iv,1)=0.E+0
19999990                                              continue
200      do 10007 i=1,n
201      do 10008 j=1,nq
202                                                   do 99989 iv=1,KDLON
203      gri(iv,j)=gri(iv,j)*exp(-tau(iv,i)/x(j))
204     &+grj(iv,i)/(alambda(iv,i)*x(j)+1.E+0)
205     &*(1.E+0-exp(-tau(iv,i)*(alambda(iv,i)+1.E+0/x(j))))
206     &+grk(iv,i)/(alambda(iv,i)*x(j)-1.E+0)
207     &*(exp(-tau(iv,i)/x(j))-exp(-tau(iv,i)*alambda(iv,i)))
208     &+sigma1(iv,i)*(1.E+0-exp(-tau(iv,i)/x(j)))
209     &+sigma2(iv,i)*(x(j)*exp(-tau(iv,i)/x(j))+tau(iv,i)-x(j))
21099989                                              continue
21110008 continue
212                                                   do 99988 iv=1,KDLON
213      fdh(iv,i+1)=0.E+0
21499988                                              continue
215      do 10009 j=1,nq
216                                                   do 99987 iv=1,KDLON
217      fdh(iv,i+1)=fdh(iv,i+1)+w(j)*x(j)*gri(iv,j)
21899987                                              continue
21910009 continue
22010007 continue
221c.......................................................................
222c on applique la condition de reflexion a sol
223c
224                                                   do 99986 iv=1,KDLON
225      fah(iv,n+1)=(1.E+0-emis(iv))*fdh(iv,n+1)+pi*emis(iv)*bsol(iv)
22699986                                              continue
227      do 10010 j=1,nq
228                                                   do 99985 iv=1,KDLON
229      gri(iv,j)=2.E+0*fah(iv,n+1)
23099985                                              continue
23110010 continue
232c.......................................................................
233c on remonte pour calculer tous les flux ascendants
234
235      do 10011 i=n,1,-1
236      do 10012 j=1,nq
237                                                   do 99984 iv=1,KDLON
238      gri(iv,j)=gri(iv,j)*exp(-tau(iv,i)/x(j))
239     &+grg(iv,i)/(alambda(iv,i)*x(j)-1.E+0)
240     &*(exp(-tau(iv,i)/x(j))-exp(-tau(iv,i)*alambda(iv,i)))
241     &+grh(iv,i)/(alambda(iv,i)*x(j)+1.E+0)
242     &*(1.E+0-exp(-tau(iv,i)*(alambda(iv,i)+1.E+0/x(j))))
243     &+alpha1(iv,i)*(1.E+0-exp(-tau(iv,i)/x(j)))
244     &+alpha2(iv,i)*(x(j)-(tau(iv,i)+x(j))*exp(-tau(iv,i)/x(j)))
24599984                                              continue
24610012 continue
247                                                   do 99983 iv=1,KDLON
248      fah(iv,i)=0.E+0
24999983                                              continue
250      do 10013 j=1,nq
251                                                   do 99982 iv=1,KDLON
252      fah(iv,i)=fah(iv,i)+w(j)*x(j)*gri(iv,j)
25399982                                              continue
25410013 continue
25510011 continue
256c.......................................................................
25710100 continue     
258c.......................................................................
259c
260c.......................................................................
261      RETURN
262      END
263
264c ***************************************************************
265
266
267      SUBROUTINE sys3v(KDLON,n,a,b,d,e,y)
268      use dimradmars_mod, only: ndlon, ndlo2, nflev
269      IMPLICIT NONE
270c.......................................................................
271c
272c  inversion d'un systeme lineaire tridiagonal
273c
274c  |  b1  d1               |   | y1 |   | e1 |
275c  |  a2  b2  d2           |   | y2 |   | e2 |
276c  |      a3  b3  d3       | * | y3 | = | e3 |
277c  |             ....      |   |    |   |    |
278c  |              an  bn   |   | yn |   | en |
279c
280c  in :   * KDLON        --> dimension de vectorisation
281c         * n            --> dimension du systeme
282c         * a,b,d,e      --> voir dessin
283c
284c  out :  * y            --> voir dessin
285c
286c.......................................................................
287c  declaration des arguments
288c
289      INTEGER KDLON,n
290      REAL a(NDLO2,n),b(NDLO2,n),d(NDLO2,n),e(NDLO2,n),y(NDLO2,n)
291c.......................................................................
292c  declaration des variables locales
293c
294      INTEGER iv,i
295      REAL as(NDLON,4*nflev),ds(NDLON,4*nflev)
296     &    ,x(NDLON,4*nflev)
297c.......................................................................
298c
299c.......................................................................
300                                                   do 99999 iv=1,KDLON
301      as(iv,n)=a(iv,n)/b(iv,n)
302      ds(iv,n)=e(iv,n)/b(iv,n)
30399999                                              continue
304      do 10001 i=n-1,1,-1
305                                                   do 99998 iv=1,KDLON
306      x(iv,i)=1.E+0/(b(iv,i)-d(iv,i)*as(iv,i+1))
307      as(iv,i)=a(iv,i)*x(iv,i)
308      ds(iv,i)=(e(iv,i)-d(iv,i)*ds(iv,i+1))*x(iv,i)
30999998                                              continue
31010001 continue
311                                                   do 99997 iv=1,KDLON
312      y(iv,1)=ds(iv,1)
31399997                                              continue
314      do 10002 i=2,n
315                                                   do 99996 iv=1,KDLON
316      y(iv,i)=ds(iv,i)-as(iv,i)*y(iv,i-1)
31799996                                              continue
31810002 continue
319c.......................................................................
320c
321c.......................................................................
322      RETURN
323      END
Note: See TracBrowser for help on using the repository browser.