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