UPP  V11.0.0
 All Data Structures Files Functions Pages
OTLIFT.f
Go to the documentation of this file.
1 
26  SUBROUTINE otlift(SLINDX)
27 
28 !
29  use vrbls3d, only: pmid, t, q
30  use vrbls2d, only: t500
31  use masks, only: lmh
32  use lookup_mod, only: thl, rdth, jtb, qs0, sqs, rdq,itb, ptbl, pl, &
33  rdp, the0, sthe, rdthe, ttbl
34  use ctlblk_mod, only: jsta, jend, im, spval, ista, iend
35  use params_mod, only: d00,h10e5, capa, elocp, eps, oneps
36  use upp_physics, only: fpvsnew
37 !
38 
39 !
40  implicit none
41 !
42 ! SET LOCAL PARAMETERS.
43  real,PARAMETER :: d8202=.820231e0 , h5e4=5.e4 , p500=50000.
44 
45 !
46 ! DECLARE VARIABLES.
47  real,intent(out) :: slindx(ista:iend,jsta:jend)
48  REAL :: tvp, esatp, qsatp
49  REAL :: tth, tp, apesp, partmp, thesp, tpsp
50  REAL :: bqs00, sqs00, bqs10, sqs10, bq, sq, tq
51  REAL :: p00, p10, p01, p11, t00, t10, t01, t11
52  REAL :: bthe00, sthe00, bthe10, sthe10, bth, sth
53  REAL :: tqq, qq, qbt, tthbt, tbt, apebt, ppq, pp
54 !
55  INTEGER :: i, j, lbtm, ittbk, iq, it, iptbk, ith, ip, iqtb
56  INTEGER :: ittb, iptb, ithtb
57 !
58 !***********************************************************************
59 ! START OTLIFT HERE
60 !
61 ! INTIALIZE LIFTED INDEX ARRAY TO ZERO.
62 !$omp parallel do private(i,j)
63  DO j=jsta,jend
64  DO i=ista,iend
65  slindx(i,j) = d00
66  ENDDO
67  ENDDO
68 !--------------FIND EXNER AT LOWEST LEVEL-------------------------------
69  DO j=jsta,jend
70  DO i=ista,iend
71  lbtm=nint(lmh(i,j))
72  IF(t(i,j,lbtm)<spval .AND. q(i,j,lbtm)<spval) THEN
73  tbt = t(i,j,lbtm)
74  qbt = q(i,j,lbtm)
75  apebt = (h10e5/pmid(i,j,lbtm))**capa
76 !--------------SCALING POTENTIAL TEMPERATURE & TABLE INDEX--------------
77  tthbt = tbt*apebt
78  tth = (tthbt-thl)*rdth
79  tqq = tth-aint(tth)
80  ittb = int(tth)+1
81 !--------------KEEPING INDICES WITHIN THE TABLE-------------------------
82  IF(ittb < 1)THEN
83  ittb = 1
84  tqq = d00
85  ENDIF
86  IF(ittb >= jtb)THEN
87  ittb = jtb-1
88  tqq = d00
89  ENDIF
90 !--------------BASE AND SCALING FACTOR FOR SPEC. HUMIDITY---------------
91  ittbk = ittb
92  bqs00=qs0(ittbk)
93  sqs00=sqs(ittbk)
94  bqs10=qs0(ittbk+1)
95  sqs10=sqs(ittbk+1)
96 !--------------SCALING SPEC. HUMIDITY & TABLE INDEX---------------------
97  bq=(bqs10-bqs00)*tqq+bqs00
98  sq=(sqs10-sqs00)*tqq+sqs00
99  tq=(qbt-bq)/sq*rdq
100  ppq = tq-aint(tq)
101  iqtb = int(tq)+1
102 !--------------KEEPING INDICES WITHIN THE TABLE-------------------------
103  IF(iqtb < 1)THEN
104  iqtb = 1
105  ppq = d00
106  ENDIF
107  IF(iqtb >= itb)THEN
108  iqtb = itb-1
109  ppq = d00
110  ENDIF
111 !--------------SATURATION PRESSURE AT FOUR SURROUNDING TABLE PTS.-------
112  iq=iqtb
113  it = ittb
114  p00=ptbl(iq,it)
115  p10=ptbl(iq+1,it)
116  p01=ptbl(iq,it+1)
117  p11=ptbl(iq+1,it+1)
118 !--------------SATURATION POINT VARIABLES AT THE BOTTOM-----------------
119  tpsp = p00+(p10-p00)*ppq+(p01-p00)*tqq &
120  +(p00-p10-p01+p11)*ppq*tqq
121  IF(tpsp <= d00) tpsp = h10e5
122  apesp = (h10e5/tpsp)**capa
123  thesp = tthbt*exp(elocp*qbt*apesp/tthbt)
124 !--------------SCALING PRESSURE & TT TABLE INDEX------------------------
125  tp=(h5e4-pl)*rdp
126  qq = tp-aint(tp)
127  iptb = int(tp)+1
128 !--------------KEEPING INDICES WITHIN THE TABLE-------------------------
129  IF(iptb < 1)THEN
130  iptb = 1
131  qq = d00
132  ENDIF
133  IF(iptb >= itb)THEN
134  iptb = itb-1
135  qq = d00
136  ENDIF
137 !--------------BASE AND SCALING FACTOR FOR THE--------------------------
138  iptbk=iptb
139  bthe00=the0(iptbk)
140  sthe00=sthe(iptbk)
141  bthe10=the0(iptbk+1)
142  sthe10=sthe(iptbk+1)
143 !--------------SCALING THE & TT TABLE INDEX-----------------------------
144  bth=(bthe10-bthe00)*qq+bthe00
145  sth=(sthe10-sthe00)*qq+sthe00
146  tth=(thesp-bth)/sth*rdthe
147  pp = tth-aint(tth)
148  ithtb = int(tth)+1
149 !--------------KEEPING INDICES WITHIN THE TABLE-------------------------
150  IF(ithtb < 1)THEN
151  ithtb = 1
152  pp = d00
153  ENDIF
154  IF(ithtb >= jtb)THEN
155  ithtb = jtb-1
156  pp = d00
157  ENDIF
158 !--------------TEMPERATURE AT FOUR SURROUNDING TT TABLE PTS.------------
159  ith=ithtb
160  ip=iptb
161  t00=ttbl(ith,ip)
162  t10=ttbl(ith+1,ip)
163  t01=ttbl(ith,ip+1)
164  t11=ttbl(ith+1,ip+1)
165 !--------------PARCEL TEMPERATURE AT 500MB----------------------------
166  IF(tpsp >= h5e4)THEN
167  partmp=(t00+(t10-t00)*pp+(t01-t00)*qq &
168  +(t00-t10-t01+t11)*pp*qq)
169  ELSE
170  partmp=tbt*apebt*d8202
171  ENDIF
172 !--------------LIFTED INDEX---------------------------------------------
173 !
174 ! GSM THE PARCEL TEMPERATURE AT 500 MB HAS BEEN COMPUTED, AND WE
175 ! FIND THE MIXING RATIO AT THAT LEVEL WHICH WILL BE THE SATURATION
176 ! VALUE SINCE WE'RE FOLLOWING A MOIST ADIABAT. NOTE THAT THE
177 ! AMBIENT 500 MB SHOULD PROBABLY BE VIRTUALIZED, BUT THE IMPACT
178 ! OF MOISTURE AT THAT LEVEL IS QUITE SMALL
179 
180  esatp=fpvsnew(partmp)
181  qsatp=eps*esatp/(p500-esatp*oneps)
182  tvp=partmp*(1+0.608*qsatp)
183  slindx(i,j)=t500(i,j)-tvp
184  ENDIF !end T(I,J,LBTM)<SPVAL
185  ENDDO
186  ENDDO
187 ! write(*,*) ' in otlift t500 partmp ',t500(1,1),partmp(1,1)
188 ! write(*,*) ' in otlift tbt ',tbt(1,1)
189 !
190  RETURN
191  END
Definition: MASKS_mod.f:1
elemental real function, public fpvsnew(t)
Definition: UPP_PHYSICS.f:345
calcape() computes CAPE/CINS and other storm related variables.
Definition: UPP_PHYSICS.f:27