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