UPP  V11.0.0
 All Data Structures Files Functions Pages
OTLFT.f
Go to the documentation of this file.
1 
28  SUBROUTINE otlft(PBND,TBND,QBND,SLINDX)
29 
30 !
31 !
32  use vrbls2d, only: t500
33  use lookup_mod, only: thl, rdth, jtb, qs0, sqs, rdq, itb, ptbl, &
34  pl, rdp, the0, sthe, rdthe, ttbl
35  use ctlblk_mod, only: jsta, jend, im, spval, ista, iend
36  use params_mod, only: d00, h10e5, capa, elocp, eps, oneps
37  use upp_physics, only: fpvsnew
38 !- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
39  implicit none
40 !
41 ! SET LOCAL PARAMETERS.
42  real,PARAMETER :: d8202=.820231e0 , h5e4=5.e4 , p500=50000.
43 
44 !
45 ! DECLARE VARIABLES.
46  real,dimension(ista:iend,jsta:jend),intent(in) :: pbnd,tbnd,qbnd
47  real,dimension(ista:iend,jsta:jend),intent(out) :: slindx
48  REAL :: tvp, esatp, qsatp
49  REAL :: bqs00, sqs00, bqs10, sqs10, p00, p10, p01, p11, bq, sq, tq
50  REAL :: bthe00, sthe00, bthe10, sthe10, bth, sth, tth
51  REAL :: t00, t10, t01, t11, tbt, qbt, apebt, tthbt, ppq, pp
52  REAL :: tqq, qq, tpsp, apesp, tthes, tp, partmp
53 !
54  INTEGER :: i, j, ittbk, iq, it, iptbk, ith, ip
55  INTEGER :: ittb, iqtb, iptb, ithtb
56 !
57 !********************************************************************
58 ! START OTLFT HERE.
59 !
60 ! ZERO LIFTED INDEX ARRAY.
61 !
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 !
69 !--------------FIND EXNER IN BOUNDARY LAYER-----------------------------
70 !
71  DO j=jsta,jend
72  DO i=ista,iend
73  tbt = tbnd(i,j)
74  qbt = qbnd(i,j)
75 !
76  if( tbt < spval ) then
77 
78  apebt = (h10e5/pbnd(i,j))**capa
79 !
80 !--------------SCALING POTENTIAL TEMPERATURE & TABLE INDEX--------------
81 !
82  tthbt = tbt*apebt
83  tth=(tthbt-thl)*rdth
84  tqq = tth-aint(tth)
85  ittb = int(tth)+1
86 !
87 !--------------KEEPING INDICES WITHIN THE TABLE-------------------------
88 !
89  IF(ittb < 1)THEN
90  ittb = 1
91  tqq = d00
92  ENDIF
93  IF(ittb >= jtb)THEN
94  ittb = jtb-1
95  tqq = d00
96  ENDIF
97 !
98 !--------------BASE AND SCALING FACTOR FOR SPEC. HUMIDITY---------------
99 !
100  ittbk = ittb
101  bqs00=qs0(ittbk)
102  sqs00=sqs(ittbk)
103  bqs10=qs0(ittbk+1)
104  sqs10=sqs(ittbk+1)
105 !
106 !--------------SCALING SPEC. HUMIDITY & TABLE INDEX---------------------
107 !
108  bq=(bqs10-bqs00)*tqq+bqs00
109  sq=(sqs10-sqs00)*tqq+sqs00
110  tq=(qbt-bq)/sq*rdq
111  ppq = tq-aint(tq)
112  iqtb = int(tq)+1
113 !
114 !--------------KEEPING INDICES WITHIN THE TABLE-------------------------
115 !
116  IF(iqtb < 1)THEN
117  iqtb = 1
118  ppq = d00
119  ENDIF
120  IF(iqtb >= itb)THEN
121  iqtb = itb-1
122  ppq = d00
123  ENDIF
124 !
125 !--------------SATURATION PRESSURE AT FOUR SURROUNDING TABLE PTS.-------
126 !
127  iq=iqtb
128  it=ittb
129  p00=ptbl(iq,it)
130  p10=ptbl(iq+1,it)
131  p01=ptbl(iq,it+1)
132  p11=ptbl(iq+1,it+1)
133 !
134 !--------------SATURATION POINT VARIABLES AT THE BOTTOM-----------------
135 !
136  tpsp = p00+(p10-p00)*ppq+(p01-p00)*tqq &
137  +(p00-p10-p01+p11)*ppq*tqq
138  IF(tpsp <= d00) tpsp = h10e5
139  apesp = (h10e5/tpsp)**capa
140  tthes = tthbt*exp(elocp*qbt*apesp/tthbt)
141 !
142 !-----------------------------------------------------------------------
143 !
144 !
145 !--------------SCALING PRESSURE & TT TABLE INDEX------------------------
146 !
147  tp = (h5e4-pl)*rdp
148  qq = tp-aint(tp)
149  iptb = int(tp)+1
150 !
151 !--------------KEEPING INDICES WITHIN THE TABLE-------------------------
152 !
153  IF(iptb < 1)THEN
154  iptb = 1
155  qq = d00
156  ENDIF
157  IF(iptb >= itb)THEN
158  iptb = itb-1
159  qq = d00
160  ENDIF
161 !
162 !--------------BASE AND SCALING FACTOR FOR THE--------------------------
163 !
164  iptbk=iptb
165  bthe00=the0(iptbk)
166  sthe00=sthe(iptbk)
167  bthe10=the0(iptbk+1)
168  sthe10=sthe(iptbk+1)
169 !
170 !--------------SCALING THE & TT TABLE INDEX-----------------------------
171 !
172  bth=(bthe10-bthe00)*qq+bthe00
173  sth=(sthe10-sthe00)*qq+sthe00
174  tth=(tthes-bth)/sth*rdthe
175  pp = tth-aint(tth)
176  ithtb = int(tth)+1
177 !
178 !--------------KEEPING INDICES WITHIN THE TABLE-------------------------
179 !
180  IF(ithtb < 1)THEN
181  ithtb = 1
182  pp = d00
183  ENDIF
184  IF(ithtb >= jtb)THEN
185  ithtb = jtb-1
186  pp = d00
187  ENDIF
188 !
189 !--------------TEMPERATURE AT FOUR SURROUNDING TT TABLE PTS.------------
190 !
191  ith=ithtb
192  ip=iptb
193  t00=ttbl(ith,ip)
194  t10=ttbl(ith+1,ip)
195  t01=ttbl(ith,ip+1)
196  t11=ttbl(ith+1,ip+1)
197 !
198 !--------------PARCEL TEMPERATURE AT 500MB----------------------------
199 !
200  IF(tpsp >= h5e4)THEN
201  partmp=(t00+(t10-t00)*pp+(t01-t00)*qq &
202  +(t00-t10-t01+t11)*pp*qq)
203  ELSE
204  partmp=tbt*apebt*d8202
205  ENDIF
206 !
207 !--------------LIFTED INDEX---------------------------------------------
208 !
209 ! GSM THE PARCEL TEMPERATURE AT 500 MB HAS BEEN COMPUTED, AND WE
210 ! FIND THE MIXING RATIO AT THAT LEVEL WHICH WILL BE THE SATURATION
211 ! VALUE SINCE WE'RE FOLLOWING A MOIST ADIABAT. NOTE THAT THE
212 ! AMBIENT 500 MB SHOULD PROBABLY BE VIRTUALIZED, BUT THE IMPACT
213 ! OF MOISTURE AT THAT LEVEL IS QUITE SMALL
214  esatp=fpvsnew(partmp)
215  qsatp=eps*esatp/(p500-esatp*oneps)
216  tvp=partmp*(1+0.608*qsatp)
217  slindx(i,j)=t500(i,j)-tvp
218 
219  else
220  slindx(i,j)=spval
221  endif
222  END DO
223  END DO
224 !
225 ! END OF ROUTINE.
226  RETURN
227  END
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