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