source: LMDZ6/branches/LMDZ_ECRad/libf/phylmdiso/cv3_enthalpmix.F90 @ 5313

Last change on this file since 5313 was 4143, checked in by dcugnet, 2 years ago
  • Some variables are renamed or replaced by direct equivalents:
    • iso_indnum -> tracers(:)%iso_iName
    • niso_possibles -> niso
    • iqiso -> iqIsoPha ; index_trac -> itZonIso
    • ok_iso_verif -> isoCheck
    • ntraceurs_zone -> nzone ; ntraciso -> ntiso
    • qperemin -> min_qparent ; masseqmin -> min_qmass ; ratiomin -> min_ratio
  • Some renamed variables are only aliased with the older name (using USE <module>, ONLY: <oldName> => <newName>) in routines where they are repeated many times.
  • Few hard-coded indexes are now computed (examples: ilic, iso, ivap, irneb, iq_vap, iq_liq, iso_H2O, iso_HDO, iso_HTO, iso_O17, iso_O18).
  • The IF(isoCheck) test is now embedded in the check_isotopes_seq and check_isotopes_loc routines (lighter calling).
  • Property svn:keywords set to Id
File size: 6.8 KB
Line 
1SUBROUTINE cv3_enthalpmix(len, nd, iflag, plim1, plim2, p, ph, &
2                       t, q, u, v, w, &
3                       wi, nk, tmix, thmix, qmix, qsmix, umix, vmix, plcl &
4#ifdef ISO
5     &                 ,xt,xtmix &
6#endif
7                       )
8#ifdef ISO
9    use infotrac_phy, ONLY: ntiso
10#endif
11  ! **************************************************************
12  ! *
13  ! CV3_ENTHALPMIX   Brassage adiabatique d'une couche d'epaisseur *
14  ! arbitraire.                                   *
15  ! *
16  ! written by   : Grandpeix Jean-Yves, 28/12/2001, 13.14.24    *
17  ! modified by :  Filiberti M-A 06/2005 vectorisation          *
18  ! **************************************************************
19
20  IMPLICIT NONE
21  ! ==============================================================
22
23  ! vertmix : determines theta, t, q, qs, u and v of the mixture generated by
24  ! adiabatic mixing of air between plim1 and plim2 with weighting w.
25  ! If plim1 and plim2 fall within the same model layer, then theta, ... v
26  ! are those of that layer.
27  ! A minimum value (dpmin) is imposed upon plim1-plim2
28
29  ! ===============================================================
30
31  include "cvthermo.h"
32  include "YOETHF.h"
33  include "YOMCST.h"
34  include "FCTTRE.h"
35!inputs:
36  INTEGER, INTENT (IN)                      :: nd, len
37  INTEGER, DIMENSION (len), INTENT (IN)     :: nk
38  REAL, DIMENSION (len), INTENT (IN)        :: plim1, plim2
39  REAL, DIMENSION (len,nd), INTENT (IN)     :: t, q
40  REAL, DIMENSION (len,nd), INTENT (IN)     :: u, v
41  REAL, DIMENSION (nd), INTENT (IN)         :: w
42  REAL, DIMENSION (len,nd), INTENT (IN)     :: p
43  REAL, DIMENSION (len,nd+1), INTENT (IN)   :: ph
44#ifdef ISO
45  REAL, DIMENSION (ntiso,len,nd), INTENT (IN)     :: xt
46#endif
47!input/output:
48  INTEGER, DIMENSION (len), INTENT (INOUT)  ::  iflag
49!outputs:
50  REAL, DIMENSION (len), INTENT (OUT)       :: tmix, thmix, qmix
51  REAL, DIMENSION (len), INTENT (OUT)       :: umix, vmix
52  REAL, DIMENSION (len), INTENT (OUT)       :: qsmix
53  REAL, DIMENSION (len), INTENT (OUT)       :: plcl
54  REAL, DIMENSION (len,nd), INTENT (OUT)    :: wi
55#ifdef ISO
56  REAL, DIMENSION (ntiso,len), INTENT (OUT)     :: xtmix
57#endif
58!internal variables :
59  INTEGER i, j
60  INTEGER niflag7
61  INTEGER, DIMENSION(len)                   :: j1, j2
62  REAL                                      :: a, b
63  REAL                                      :: cpn
64  REAL                                      :: x, y, p0, p0m1, zdelta, zcor
65  REAL, SAVE                                :: dpmin=1.
66!$OMP THREADPRIVATE(dpmin)
67  REAL, DIMENSION(len)                      :: plim2p  ! = min(plim2(:),plim1(:)-dpmin)
68  REAL, DIMENSION(len)                      :: akm     ! mixture enthalpy
69  REAL, DIMENSION(len)                      :: dpw, coef
70  REAL, DIMENSION(len)                      :: rdcp, a2, b2, pnk
71  REAL, DIMENSION(len)                      :: rh, chi
72  REAL, DIMENSION(len)                      :: eqwght
73  REAL, DIMENSION(len,nd)                   :: p1, p2
74#ifdef ISO
75  INTEGER ixt
76#endif
77
78
79!!  print *,' ->cv3_vertmix, plim1,plim2 ', plim1,plim2   !jyg
80  plim2p(:) = min(plim2(:),plim1(:)-dpmin)
81  j1(:)=nd
82  j2(:) = 0
83  DO j = 1, nd
84    DO i = 1, len
85      IF (plim1(i)<=ph(i,j)) j1(i) = j
86!!!      IF (plim2p(i)>=ph(i,j+1) .AND. plim2p(i)<ph(i,j)) j2(i) = j
87      IF (plim2p(i)< ph(i,j)) j2(i) = j
88    END DO
89  END DO
90
91  DO j = 1, nd
92    DO i = 1, len
93      wi(i, j) = 0.
94    END DO
95  END DO
96  DO i = 1, len
97    akm(i) = 0.
98    qmix(i) = 0.
99    umix(i) = 0.
100    vmix(i) = 0.
101    dpw(i) = 0.
102    a2(i) = 0.0
103    b2(i) = 0.
104    pnk(i) = p(i, nk(i))
105#ifdef ISO
106    xtmix(:,i) = 0.
107#endif
108  END DO
109  eqwght(:) = 0.
110
111  p0 = 1000.
112  p0m1 = 1./p0
113
114  DO i = 1, len
115    IF (j2(i) < j1(i)) THEN
116      coef(i) = 1.
117      eqwght(i) = 1.
118    ELSE
119      coef(i) = 1./(plim1(i)-plim2p(i))
120    ENDIF
121  END DO
122
123!!  print *,'cv3_vertmix, j1,j2,coef ', j1,j2,coef  !jyg
124
125  DO j = 1, nd
126    DO i = 1, len
127      IF (j>=j1(i) .AND. j<=j2(i)) THEN
128        p1(i, j) = min(ph(i,j), plim1(i))
129        p2(i, j) = max(ph(i,j+1), plim2p(i))
130        ! CRtest:couplage thermiques: deja normalise
131        ! wi(i,j) = w(j)
132        ! print*,'wi',wi(i,j)
133        wi(i, j) = w(j)*(p1(i,j)-p2(i,j))*coef(i)+eqwght(i)
134        dpw(i) = dpw(i) + wi(i, j)
135
136!!  print *,'cv3_vertmix, j, wi(1,j),dpw ', j, wi(1,j),dpw  !jyg
137
138      END IF
139    END DO
140  END DO
141
142  ! CR:print
143  ! do i=1,len
144  ! print*,'plim',plim1(i),plim2p(i)
145  ! enddo
146  DO j = 1, nd
147    DO i = 1, len
148      IF (j>=j1(i) .AND. j<=j2(i)) THEN
149        wi(i, j) = wi(i, j)/dpw(i)
150        akm(i) = akm(i) + (cpd*(1.-q(i,j))+q(i,j)*cpv)*t(i, j)*wi(i, j)
151        qmix(i) = qmix(i) + q(i, j)*wi(i, j)
152        umix(i) = umix(i) + u(i, j)*wi(i, j)
153        vmix(i) = vmix(i) + v(i, j)*wi(i, j)
154#ifdef ISO
155        do ixt=1,ntiso
156          xtmix(ixt,i) = xtmix(ixt,i) +  xt(ixt,i, j)*wi(i, j)
157        enddo
158#endif
159      END IF
160    END DO
161  END DO
162
163  DO i = 1, len
164    rdcp(i) = (rrd*(1.-qmix(i))+qmix(i)*rrv)/(cpd*(1.-qmix(i))+qmix(i)*cpv)
165  END DO
166
167
168!!  print *,'cv3_vertmix, rdcp ', rdcp  !jyg
169
170
171
172  DO j = 1, nd
173    DO i = 1, len
174      IF (j>=j1(i) .AND. j<=j2(i)) THEN
175        ! c            x=(.5*(p1(i,j)+p2(i,j))*p0m1)**rdcp(i)
176        y = (.5*(p1(i,j)+p2(i,j))/pnk(i))**rdcp(i)
177        ! c            a2(i)=a2(i)+(cpd*(1.-qmix(i))+qmix(i)*cpv)*x*wi(i,j)
178        b2(i) = b2(i) + (cpd*(1.-qmix(i))+qmix(i)*cpv)*y*wi(i, j)
179      END IF
180    END DO
181  END DO
182
183  DO i = 1, len
184    tmix(i) = akm(i)/b2(i)
185    thmix(i) = tmix(i)*(p0/pnk(i))**rdcp(i)
186    ! print*,'thmix akm',akm(i),b2(i)
187    ! print*,'thmix t',tmix(i),p0
188    ! print*,'thmix p',pnk(i),rdcp(i)
189    ! print*,'thmix',thmix(i)
190    ! c         thmix(i) = akm(i)/a2(i)
191    ! c         tmix(i)= thmix(i)*(pnk(i)*p0m1)**rdcp(i)
192    zdelta = max(0., sign(1.,rtt-tmix(i)))
193    qsmix(i) = r2es*foeew(tmix(i), zdelta)/(pnk(i)*100.)
194    qsmix(i) = min(0.5, qsmix(i))
195    zcor = 1./(1.-retv*qsmix(i))
196    qsmix(i) = qsmix(i)*zcor
197  END DO
198
199  ! -------------------------------------------------------------------
200  ! --- Calculate lifted condensation level of air at parcel origin level
201  ! --- (Within 0.2% of formula of Bolton, MON. WEA. REV.,1980)
202  ! -------------------------------------------------------------------
203
204  a = 1669.0 ! convect3
205  b = 122.0 ! convect3
206
207
208  niflag7 = 0
209  DO i = 1, len
210
211    IF (iflag(i)/=7) THEN ! modif sb Jun7th 2002
212
213      rh(i) = qmix(i)/qsmix(i)
214      chi(i) = tmix(i)/(a-b*rh(i)-tmix(i)) ! convect3
215      ! ATTENTION, la LIGNE DESSOUS A ETE RAJOUTEE ARBITRAIREMENT ET
216      ! MASQUE UN PB POTENTIEL
217      chi(i) = max(chi(i), 0.)
218      rh(i) = max(rh(i), 0.)
219      plcl(i) = pnk(i)*(rh(i)**chi(i))
220      IF (((plcl(i)<200.0) .OR. (plcl(i)>=2000.0)) .AND. (iflag(i)==0)) &
221          iflag(i) = 8
222
223    ELSE
224
225      niflag7 = niflag7 + 1
226      plcl(i) = plim2p(i)
227
228    END IF ! iflag=7
229
230    ! print*,'NIFLAG7  =',niflag7
231
232  END DO
233
234!!  print *,' cv3_vertmix->'  !jyg
235
236
237  RETURN
238END SUBROUTINE cv3_enthalpmix
239
Note: See TracBrowser for help on using the repository browser.