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

Last change on this file since 317 was 38, checked in by emillour, 14 years ago

Ajout du modè Martien (mon LMDZ.MARS.BETA, du 28/01/2011) dans le rértoire mars, pour pouvoir suivre plus facilement les modifs.
EM

File size: 13.9 KB
Line 
1      SUBROUTINE flusv(KDLON,nsf,n,omega,g,tau,emis,bh,bsol,fah,fdh)
2      IMPLICIT NONE
3c.......................................................................
4c
5c  calcul des flux ascendant et descendant aux interfaces entre n couches
6c  * dans l'infrarouge
7c  * B est une fonction lineaire de $\tau$ a l'interieur de chaque couche
8c  * le B du sol peut etre different de celui qui correspond au profil
9c    de la n-ieme couche
10c  * l'hypothese est une hypothese a deux flux isotropes sur chaque
11c    hemisphere ("hemispheric constant") + "source function technique"
12c    (voir Toon et al. 1988)
13c  * le flux descendant en haut de l'atmosphere est nul
14c  * les couches sont numerotees du haut de l'atmosphere vers le sol
15c
16c  in :   * KDLON      ---> dimension de vectorisation
17c         * nsf        ---> nsf=0 ==> "hemispheric constant"
18c                           nsf>0 ==> "hemispheric constant" + "source function"
19c         * n          ---> nombre de couches
20c         * omega(i)   ---> single scattering albedo pour la i-eme couche
21c         * g(i)       ---> asymmetry parameter pour la i-eme couche
22c         * tau(i)     ---> epaisseur optique de la i-eme couche
23c         * emis       ---> emissivite du sol
24c         * bh(i)      ---> luminance du corps noir en haut de la i-eme
25c                           couche, bh(n+1) pour la valeur au sol qui
26c                           correspond au profil de la n-ieme couche
27c         * bsol       ---> luminance du corps noir au sol
28c
29c  out :  * fah(i)     ---> flux ascendant en haut de la i-eme couche,
30c                           fah(n+1) pour le sol
31c         * fdh(i)     ---> flux descendant en haut de la i-eme couche,
32c                           fdh(n+1) pour le sol
33c
34c.......................................................................
35c  include des dimensions locales
36c
37#include "dimensions.h"
38#include "dimphys.h"
39#include "dimradmars.h"
40c.......................................................................
41c  declaration des arguments
42c
43      INTEGER KDLON,nsf,n
44      REAL omega(NDLO2,n),g(NDLO2,n),tau(NDLO2,n),emis(NDLO2)
45     &,bh(NDLO2,n+1),bsol(NDLO2),fah(NDLO2,n+1),fdh(NDLO2,n+1)
46c.......................................................................
47c  declaration des variables locales
48c
49      REAL pi
50      PARAMETER (pi=3.141592653589793E+0)
51      INTEGER iv,i,j
52      REAL beta,gama1,gama2,amu1,grgama,b0,b1
53      REAL a(NDLON,4*nlayermx),b(NDLON,4*nlayermx)
54     &    ,d(NDLON,4*nlayermx),e(NDLON,4*nlayermx)
55     &    ,y(NDLON,4*nlayermx)
56     &    ,alambda(NDLON,2*nlayermx)
57     &    ,e1(NDLON,2*nlayermx),e2(NDLON,2*nlayermx)
58     &    ,e3(NDLON,2*nlayermx),e4(NDLON,2*nlayermx)
59     &    ,cah(NDLON,2*nlayermx),cab(NDLON,2*nlayermx)
60     &    ,cdh(NDLON,2*nlayermx),cdb(NDLON,2*nlayermx)
61      REAL grg(NDLON,2*nlayermx),grh(NDLON,2*nlayermx)
62     &    ,grj(NDLON,2*nlayermx),grk(NDLON,2*nlayermx)
63     &    ,alpha1(NDLON,2*nlayermx),alpha2(NDLON,2*nlayermx)
64     &    ,sigma1(NDLON,2*nlayermx),sigma2(NDLON,2*nlayermx)
65      INTEGER nq
66      PARAMETER (nq=8)
67      REAL x(nq),w(nq),gri(NDLON,nq)
68      DATA x/1.9855071751231860E-2 , 0.1016667612931866E+0
69     &     , 0.2372337950418355E+0 , 0.4082826787521751E+0
70     &     , 0.5917173212478250E+0 , 0.7627662049581645E+0
71     &     , 0.8983332387068134E+0 , 0.9801449282487682E+0/
72      DATA w/5.0614268145185310E-2 , 0.1111905172266872E+0
73     &     , 0.1568533229389437E+0 , 0.1813418916891810E+0
74     &     , 0.1813418916891810E+0 , 0.1568533229389437E+0
75     &     , 0.1111905172266872E+0 , 5.0614268145185310E-2/
76c.......................................................................
77c
78c.......................................................................
79      do 10001 i=1,n
80                                                   do 99999 iv=1,KDLON
81      beta=(1.E+0-g(iv,i))/2.E+0
82      gama1=2.E+0*(1.E+0-omega(iv,i)*(1.E+0-beta))
83      gama2=2.E+0*omega(iv,i)*beta
84      amu1=5.E-1
85      alambda(iv,i)=sqrt(gama1**2-gama2**2)
86      grgama=(gama1-alambda(iv,i))/gama2
87c
88c ici une petite bidouille : si l'epaisseur optique d'une couche
89c est trop faible, $dB \over d\tau$ devient tres grand et le schema
90c s'ecroule. dans ces cas, on fait l'hypothese de couche isotherme.
91c
92      if (tau(iv,i).gt.1.E-3) then
93      b0=bh(iv,i)
94      b1=(bh(iv,i+1)-b0)/tau(iv,i)
95      else
96      b0=(bh(iv,i)+bh(iv,i+1))/2.E+0
97      b1=0.E+0
98      endif
99c
100      e1(iv,i)=1.E+0+grgama*exp(-alambda(iv,i)*tau(iv,i))
101      e2(iv,i)=1.E+0-grgama*exp(-alambda(iv,i)*tau(iv,i))
102      e3(iv,i)=grgama+exp(-alambda(iv,i)*tau(iv,i))
103      e4(iv,i)=grgama-exp(-alambda(iv,i)*tau(iv,i))
104      cah(iv,i)=2.E+0*pi*amu1*(b0+b1/(gama1+gama2))
105      cab(iv,i)=2.E+0*pi*amu1*(b0+b1*(tau(iv,i)+1.E+0/(gama1+gama2)))
106      cdh(iv,i)=2.E+0*pi*amu1*(b0-b1/(gama1+gama2))
107      cdb(iv,i)=2.E+0*pi*amu1*(b0+b1*(tau(iv,i)-1.E+0/(gama1+gama2))) 
108c
109      grg(iv,i)=(1.E+0/amu1-alambda(iv,i))
110      grh(iv,i)=grgama*(alambda(iv,i)+1.E+0/amu1)
111      grj(iv,i)=grh(iv,i)
112      grk(iv,i)=grg(iv,i)
113      alpha1(iv,i)=2.E+0*pi*(b0+b1*(1.E+0/(gama1+gama2)-amu1))
114      alpha2(iv,i)=2.E+0*pi*b1
115      sigma1(iv,i)=2.E+0*pi*(b0-b1*(1.E+0/(gama1+gama2)-amu1))
116      sigma2(iv,i)=alpha2(iv,i)
117c
11899999                                              continue
11910001 continue
120c.......................................................................
121                                                   do 99998 iv=1,KDLON
122      a(iv,1)=0.E+0
123      b(iv,1)=e1(iv,1)
124      d(iv,1)=-e2(iv,1)
125      e(iv,1)=-cdh(iv,1)
12699998                                              continue
127c
128      do 10002 i=1,n-1
129      j=2*i+1
130                                                   do 99997 iv=1,KDLON
131      a(iv,j)=e2(iv,i)*e3(iv,i)-e4(iv,i)*e1(iv,i)
132      b(iv,j)=e1(iv,i)*e1(iv,i+1)-e3(iv,i)*e3(iv,i+1)
133      d(iv,j)=e3(iv,i)*e4(iv,i+1)-e1(iv,i)*e2(iv,i+1)
134      e(iv,j)=e3(iv,i)*(cah(iv,i+1)-cab(iv,i))
135     &       +e1(iv,i)*(cdb(iv,i)-cdh(iv,i+1))
13699997                                              continue
13710002 continue
138c
139      do 10003 i=1,n-1
140      j=2*i
141                                                   do 99996 iv=1,KDLON
142      a(iv,j)=e2(iv,i+1)*e1(iv,i)-e3(iv,i)*e4(iv,i+1)
143      b(iv,j)=e2(iv,i)*e2(iv,i+1)-e4(iv,i)*e4(iv,i+1)
144      d(iv,j)=e1(iv,i+1)*e4(iv,i+1)-e2(iv,i+1)*e3(iv,i+1)
145      e(iv,j)=e2(iv,i+1)*(cah(iv,i+1)-cab(iv,i))
146     &       +e4(iv,i+1)*(cdb(iv,i)-cdh(iv,i+1))
14799996                                              continue
14810003 continue
149c     
150      j=2*n
151                                                   do 99995 iv=1,KDLON
152      a(iv,j)=e1(iv,n)-(1.E+0-emis(iv))*e3(iv,n)
153      b(iv,j)=e2(iv,n)-(1.E+0-emis(iv))*e4(iv,n)
154      d(iv,j)=0.E+0
155      e(iv,j)=emis(iv)*pi*bsol(iv)-cab(iv,n)+(1.E+0-emis(iv))*cdb(iv,n)
15699995                                              continue
157c.......................................................................
158      call sys3v(KDLON,2*n,a,b,d,e,y)
159c.......................................................................
160      do 10004 i=1,n
161                                                   do 99994 iv=1,KDLON
162      grg(iv,i)=grg(iv,i)*(y(iv,2*i-1)+y(iv,2*i))
163      grh(iv,i)=grh(iv,i)*(y(iv,2*i-1)-y(iv,2*i))
164      grj(iv,i)=grj(iv,i)*(y(iv,2*i-1)+y(iv,2*i))
165      grk(iv,i)=grk(iv,i)*(y(iv,2*i-1)-y(iv,2*i))
16699994                                              continue
16710004 continue
168c.......................................................................
169c les valeurs de flux "hemispheric constant"
170c
171      IF (nsf.eq.0) THEN
172      do 10005 i=1,n
173                                                   do 99993 iv=1,KDLON
174      fah(iv,i)=e3(iv,i)*y(iv,2*i-1)-e4(iv,i)*y(iv,2*i)+cah(iv,i)
175      fdh(iv,i)=e1(iv,i)*y(iv,2*i-1)-e2(iv,i)*y(iv,2*i)+cdh(iv,i)
17699993                                              continue
17710005 continue
178                                                   do 99992 iv=1,KDLON
179      fah(iv,n+1)=e1(iv,n)*y(iv,2*n-1)+e2(iv,n)*y(iv,2*n)+cab(iv,n)
180      fdh(iv,n+1)=e3(iv,n)*y(iv,2*n-1)+e4(iv,n)*y(iv,2*n)+cdb(iv,n)
18199992                                              continue
182      GOTO 10100
183      ENDIF
184c.......................................................................
185c passage a "source function"
186c
187c on applique une quadrature de dimension nq
188c x est le vecteur des \mu de la quadrature
189c w est le vecteur des poids correspondants
190c nq est fixe en parameter
191c x et w sont fixes en data
192c
193c.......................................................................
194c on part d'en haut et on descent selon les nq angles pour calculer
195c tous les flux descendants
196c
197      do 10006 j=1,nq
198                                                   do 99991 iv=1,KDLON
199      gri(iv,j)=0.E+0
20099991                                              continue
20110006 continue
202                                                   do 99990 iv=1,KDLON
203      fdh(iv,1)=0.E+0
20499990                                              continue
205      do 10007 i=1,n
206      do 10008 j=1,nq
207                                                   do 99989 iv=1,KDLON
208      gri(iv,j)=gri(iv,j)*exp(-tau(iv,i)/x(j))
209     &+grj(iv,i)/(alambda(iv,i)*x(j)+1.E+0)
210     &*(1.E+0-exp(-tau(iv,i)*(alambda(iv,i)+1.E+0/x(j))))
211     &+grk(iv,i)/(alambda(iv,i)*x(j)-1.E+0)
212     &*(exp(-tau(iv,i)/x(j))-exp(-tau(iv,i)*alambda(iv,i)))
213     &+sigma1(iv,i)*(1.E+0-exp(-tau(iv,i)/x(j)))
214     &+sigma2(iv,i)*(x(j)*exp(-tau(iv,i)/x(j))+tau(iv,i)-x(j))
21599989                                              continue
21610008 continue
217                                                   do 99988 iv=1,KDLON
218      fdh(iv,i+1)=0.E+0
21999988                                              continue
220      do 10009 j=1,nq
221                                                   do 99987 iv=1,KDLON
222      fdh(iv,i+1)=fdh(iv,i+1)+w(j)*x(j)*gri(iv,j)
22399987                                              continue
22410009 continue
22510007 continue
226c.......................................................................
227c on applique la condition de reflexion a sol
228c
229                                                   do 99986 iv=1,KDLON
230      fah(iv,n+1)=(1.E+0-emis(iv))*fdh(iv,n+1)+pi*emis(iv)*bsol(iv)
23199986                                              continue
232      do 10010 j=1,nq
233                                                   do 99985 iv=1,KDLON
234      gri(iv,j)=2.E+0*fah(iv,n+1)
23599985                                              continue
23610010 continue
237c.......................................................................
238c on remonte pour calculer tous les flux ascendants
239
240      do 10011 i=n,1,-1
241      do 10012 j=1,nq
242                                                   do 99984 iv=1,KDLON
243      gri(iv,j)=gri(iv,j)*exp(-tau(iv,i)/x(j))
244     &+grg(iv,i)/(alambda(iv,i)*x(j)-1.E+0)
245     &*(exp(-tau(iv,i)/x(j))-exp(-tau(iv,i)*alambda(iv,i)))
246     &+grh(iv,i)/(alambda(iv,i)*x(j)+1.E+0)
247     &*(1.E+0-exp(-tau(iv,i)*(alambda(iv,i)+1.E+0/x(j))))
248     &+alpha1(iv,i)*(1.E+0-exp(-tau(iv,i)/x(j)))
249     &+alpha2(iv,i)*(x(j)-(tau(iv,i)+x(j))*exp(-tau(iv,i)/x(j)))
25099984                                              continue
25110012 continue
252                                                   do 99983 iv=1,KDLON
253      fah(iv,i)=0.E+0
25499983                                              continue
255      do 10013 j=1,nq
256                                                   do 99982 iv=1,KDLON
257      fah(iv,i)=fah(iv,i)+w(j)*x(j)*gri(iv,j)
25899982                                              continue
25910013 continue
26010011 continue
261c.......................................................................
26210100 continue     
263c.......................................................................
264c
265c.......................................................................
266      RETURN
267      END
268
269c ***************************************************************
270
271
272      SUBROUTINE sys3v(KDLON,n,a,b,d,e,y)
273      IMPLICIT NONE
274c.......................................................................
275c
276c  inversion d'un systeme lineaire tridiagonal
277c
278c  |  b1  d1               |   | y1 |   | e1 |
279c  |  a2  b2  d2           |   | y2 |   | e2 |
280c  |      a3  b3  d3       | * | y3 | = | e3 |
281c  |             ....      |   |    |   |    |
282c  |              an  bn   |   | yn |   | en |
283c
284c  in :   * KDLON        --> dimension de vectorisation
285c         * n            --> dimension du systeme
286c         * a,b,d,e      --> voir dessin
287c
288c  out :  * y            --> voir dessin
289c
290c.......................................................................
291c  include des dimensions locales
292c
293#include "dimensions.h"
294#include "dimphys.h"
295#include "dimradmars.h"
296c.......................................................................
297c  declaration des arguments
298c
299      INTEGER KDLON,n
300      REAL a(NDLO2,n),b(NDLO2,n),d(NDLO2,n),e(NDLO2,n),y(NDLO2,n)
301c.......................................................................
302c  declaration des variables locales
303c
304      INTEGER iv,i
305      REAL as(NDLON,4*nlayermx),ds(NDLON,4*nlayermx)
306     &    ,x(NDLON,4*nlayermx)
307c.......................................................................
308c
309c.......................................................................
310                                                   do 99999 iv=1,KDLON
311      as(iv,n)=a(iv,n)/b(iv,n)
312      ds(iv,n)=e(iv,n)/b(iv,n)
31399999                                              continue
314      do 10001 i=n-1,1,-1
315                                                   do 99998 iv=1,KDLON
316      x(iv,i)=1.E+0/(b(iv,i)-d(iv,i)*as(iv,i+1))
317      as(iv,i)=a(iv,i)*x(iv,i)
318      ds(iv,i)=(e(iv,i)-d(iv,i)*ds(iv,i+1))*x(iv,i)
31999998                                              continue
32010001 continue
321                                                   do 99997 iv=1,KDLON
322      y(iv,1)=ds(iv,1)
32399997                                              continue
324      do 10002 i=2,n
325                                                   do 99996 iv=1,KDLON
326      y(iv,i)=ds(iv,i)-as(iv,i)*y(iv,i-1)
32799996                                              continue
32810002 continue
329c.......................................................................
330c
331c.......................................................................
332      RETURN
333      END
Note: See TracBrowser for help on using the repository browser.