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