source: LMDZ4/branches/LMDZ4V5.0-LF/libf/phylmd/cltrac.F90 @ 5427

Last change on this file since 5427 was 1279, checked in by Laurent Fairhead, 15 years ago

Merged LMDZ4-dev branch changes r1241:1278 into the trunk
Running trunk and LMDZ4-dev in LMDZOR configuration on local
machine (sequential) and SX8 (4-proc) yields identical results
(restart and restartphy are identical binarily)
Log history from r1241 to r1278 is available by switching to
source:LMDZ4/branches/LMDZ4-dev-20091210

  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 4.3 KB
Line 
1!
2! $Id $
3!
4SUBROUTINE cltrac(dtime,coef,t,tr,flux,paprs,pplay,delp,d_tr)
5  USE dimphy
6  IMPLICIT NONE
7!======================================================================
8! Auteur(s): O. Boucher (LOA/LMD) date: 19961127
9!            inspire de clvent
10! Objet: diffusion verticale de traceurs avec flux fixe a la surface
11!        ou/et flux du type c-drag
12!
13! Arguments:
14!-----------
15! dtime....input-R- intervalle du temps (en secondes)
16! coef.....input-R- le coefficient d'echange (m**2/s) l>1
17! t........input-R- temperature (K)
18! tr.......input-R- la q. de traceurs
19! flux.....input-R- le flux de traceurs a la surface
20! paprs....input-R- pression a inter-couche (Pa)
21! pplay....input-R- pression au milieu de couche (Pa)
22! delp.....input-R- epaisseur de couche (Pa)
23! cdrag....input-R- cdrag pour le flux de surface (non active)
24! tr0......input-R- traceurs a la surface ou dans l'ocean (non active)
25! d_tr.....output-R- le changement de tr
26! flux_tr..output-R- flux de tr
27!======================================================================
28  include "YOMCST.h"
29!
30! Entree
31!
32  REAL,INTENT(IN)                        :: dtime
33  REAL,DIMENSION(klon,klev),INTENT(IN)   :: coef
34  REAL,DIMENSION(klon,klev),INTENT(IN)   :: t, tr
35  REAL,DIMENSION(klon),INTENT(IN)        :: flux !(at/s/m2)
36  REAL,DIMENSION(klon,klev+1),INTENT(IN) :: paprs
37  REAL,DIMENSION(klon,klev),INTENT(IN)   :: pplay, delp
38!
39! Sorties
40!
41  REAL ,DIMENSION(klon,klev),INTENT(OUT) :: d_tr
42!  REAL ,DIMENSION(klon,klev),INTENT(OUT) :: flux_tr
43!
44! Local
45!
46  INTEGER                   :: i, k
47  REAL,DIMENSION(klon)      :: cdrag, tr0
48  REAL,DIMENSION(klon,klev) :: zx_ctr
49  REAL,DIMENSION(klon,klev) :: zx_dtr
50  REAL,DIMENSION(klon)      :: zx_buf
51  REAL,DIMENSION(klon,klev) :: zx_coef
52  REAL,DIMENSION(klon,klev) :: local_tr
53  REAL,DIMENSION(klon)      :: zx_alf1,zx_alf2,zx_flux
54
55!======================================================================
56
57  DO k = 1, klev
58     DO i = 1, klon
59        local_tr(i,k) = tr(i,k)
60     ENDDO
61  ENDDO
62
63!======================================================================
64
65  DO i = 1, klon
66     zx_alf1(i) = (paprs(i,1)-pplay(i,2))/(pplay(i,1)-pplay(i,2))
67     zx_alf2(i) = 1.0 - zx_alf1(i)
68     zx_flux(i) =  -flux(i)*dtime*RG
69! Pour le moment le flux est prescrit cdrag et zx_coef(1) vaut 0
70     cdrag(i) = 0.0
71     tr0(i) = 0.0
72     zx_coef(i,1) = cdrag(i)*dtime*RG
73     zx_ctr(i,1)=0.
74     zx_dtr(i,1)=0.
75  ENDDO
76
77!======================================================================
78
79  DO k = 2, klev
80     DO i = 1, klon
81        zx_coef(i,k) = coef(i,k)*RG/(pplay(i,k-1)-pplay(i,k))   &
82             *(paprs(i,k)*2/(t(i,k)+t(i,k-1))/RD)**2
83        zx_coef(i,k) = zx_coef(i,k)*dtime*RG 
84     ENDDO
85  ENDDO
86
87!======================================================================
88
89  DO i = 1, klon
90     zx_buf(i) = delp(i,1) + zx_coef(i,1)*zx_alf1(i) + zx_coef(i,2)
91     !
92     zx_ctr(i,2) = (local_tr(i,1)*delp(i,1)+                  &
93          zx_coef(i,1)*tr0(i)-zx_flux(i))/zx_buf(i)
94     !
95     zx_dtr(i,2) = (zx_coef(i,2)-zx_alf2(i)*zx_coef(i,1)) /   &
96          zx_buf(i)
97  ENDDO
98
99  DO k = 3, klev
100     DO i = 1, klon
101        zx_buf(i) = delp(i,k-1) + zx_coef(i,k)      &
102             + zx_coef(i,k-1)*(1.-zx_dtr(i,k-1))
103        zx_ctr(i,k) = (local_tr(i,k-1)*delp(i,k-1)  &
104             +zx_coef(i,k-1)*zx_ctr(i,k-1) )/zx_buf(i)
105        zx_dtr(i,k) = zx_coef(i,k)/zx_buf(i)
106     ENDDO
107  ENDDO
108
109  DO i = 1, klon
110     local_tr(i,klev) = ( local_tr(i,klev)*delp(i,klev) &
111          +zx_coef(i,klev)*zx_ctr(i,klev) )             &
112          / ( delp(i,klev) + zx_coef(i,klev)            &
113          -zx_coef(i,klev)*zx_dtr(i,klev) )
114  ENDDO
115
116  DO k = klev-1, 1, -1
117     DO i = 1, klon
118        local_tr(i,k) = zx_ctr(i,k+1) + zx_dtr(i,k+1)*local_tr(i,k+1)
119     ENDDO
120  ENDDO
121
122!======================================================================
123!== flux_tr est le flux de traceur (positif vers bas)
124!      DO i = 1, klon
125!         flux_tr(i,1) = zx_coef(i,1)/(RG*dtime)
126!      ENDDO
127!      DO k = 2, klev
128!      DO i = 1, klon
129!         flux_tr(i,k) = zx_coef(i,k)/(RG*dtime)
130!     .               * (local_tr(i,k)-local_tr(i,k-1))
131!      ENDDO
132!      ENDDO
133!======================================================================
134  DO k = 1, klev
135     DO i = 1, klon
136        d_tr(i,k) = local_tr(i,k) - tr(i,k)
137     ENDDO
138  ENDDO
139 
140END SUBROUTINE cltrac
Note: See TracBrowser for help on using the repository browser.