vcoord_gen 1.14.0
Loading...
Searching...
No Matches
vcoord_gen.f90
Go to the documentation of this file.
1
4
151subroutine vcoord_gen(levs,lupp,pbot,psig,ppre,pupp,ptop,&
152 dpbot,dpsig,dppre,dpupp,dptop,pmin,ak,bk)
153 implicit none
154 integer,intent(in):: levs,lupp
155 real,intent(in):: pbot,psig,ppre,pupp,ptop
156 real,intent(in):: dpbot,dpsig,dppre,dpupp,dptop
157 real,intent(out):: pmin,ak(0:levs),bk(0:levs)
158 integer,parameter:: lo=100,li=10 ! outer and inner N-R iterations
159 real pdif ! thickness from pbot to pupp
160 real delb ! delta sig at bot
161 real delt ! delta sig at top
162 real sig1 ! crossing parameter 1
163 real sig2 ! crossing parameter 2
164 real c1 ! proportionality constant
165 real c2 ! integration constant
166 real sig ! sig variable
167 real dsig ! delta sig variable
168 real delbio0 ! initial guess at delta sig at bot
169 real deltio0 ! initial guess at delta sig at top
170 real delbio ! guess at delta sig at bot
171 real deltio ! guess at delta sig at top
172 real c1sig1 ! d(c1)/d(sig1)
173 real c1sig2 ! d(c1)/d(sig2)
174 real p(2) ! rhs in N-R iteration
175 real fjac(2,2) ! lhs in N-R iteration
176 integer indx(2) ! permutations in N-R iteration
177 real ppred ! pressure at T1-T2 boundary
178 real spre ! sig at P-T1 boundary
179 real spred ! sig at T1-T2 boundary
180 real rkpre ! level at P-T1 boundary
181 real rkpred ! level at T1-T2 boundary
182 real pkpre ! dp/dk at P-T1 boundary
183 real pkkpre ! d2p/dk2 at P-T1 boundary
184 real psigd ! pressure at T2-T3 boundary
185 real ssig ! sig at T3-S boundary
186 real ssigd ! sig at T2-T3 boundary
187 real rksig ! level at T3-S boundary
188 real rksigd ! level at T2-T3 boundary
189 real pksig ! dp/dk at T3-S boundary
190 real pkksig ! d2p/dk2 at T3-S boundary
191 real p2sig ! pressure at T3-S boundary for pmin surface pressure
192 real p2sigd ! pressure at T2-T3 boundary for pmin surface pressure
193 real p2pred ! pressure at T1-T2 boundary for pmin surface pressure
194 real x2(9) ! rhs in linear solver
195 real a2(9,9) ! lhs in linear solver
196 integer indx2(9) ! permutations in linear solver
197 real pkupp ! dp/dk at U-P boundary
198 real pkkupp ! d2p/dk2 at U-P boundary
199 real x3(3) ! rhs in linear solver
200 real a3(3,3) ! lhs in linear solver
201 integer indx3(3) ! permutations in linear solver
202 real p1 ! pressure variable for pbot surface pressure
203 real p2 ! pressure variable for pmin surface pressure
204 real d ! determinant permutation
205 integer io,ii,k
206! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
207! STEP 1.
208 pdif=pbot-pupp
209 delb=dpbot/pdif
210 delt=dpupp/pdif
211 sig1=1+delb
212 sig2=-delt
213 c1=log(-sig2*(1-sig1)/sig1/(sig2-1))/lupp
214 c2=log((sig2-1)/(1-sig1))
215 sig=1
216 dsig=(sig2-sig)*(sig-sig1)*c1/(sig1-sig2)
217 delbio0=-dsig
218 sig=0
219 dsig=(sig2-sig)*(sig-sig1)*c1/(sig1-sig2)
220 deltio0=-dsig
221! Newton-Raphson iterations
222 do io=1,lo
223 delbio=delbio0+(delb-delbio0)*io/lo
224 deltio=deltio0+(delt-deltio0)*io/lo
225 do ii=1,li
226 c1sig1=-1/(sig1*(1-sig1)*lupp)
227 c1sig2=-1/(sig2*(sig2-1)*lupp)
228 sig=1
229 dsig=(sig2-sig)*(sig-sig1)*c1/(sig1-sig2)
230 p(1)=-delbio-dsig
231 fjac(1,1)=((-c1*(sig+sig2)+(sig-sig1)*c1sig1*(sig1+sig2)) &
232 *(sig2-sig)/(sig1+sig2)**2)
233 fjac(1,2)=((c1*(sig+sig1)+(sig2-sig)*c1sig2*(sig1+sig2)) &
234 *(sig-sig1)/(sig1+sig2)**2)
235 sig=0
236 dsig=(sig2-sig)*(sig-sig1)*c1/(sig1-sig2)
237 p(2)=-deltio-dsig
238 fjac(2,1)=((-c1*(sig+sig2)+(sig-sig1)*c1sig1*(sig1+sig2)) &
239 *(sig2-sig)/(sig1+sig2)**2)
240 fjac(2,2)=((c1*(sig+sig1)+(sig2-sig)*c1sig2*(sig1+sig2)) &
241 *(sig-sig1)/(sig1+sig2)**2)
242 call ludcmp(fjac,2,2,indx,d)
243 call lubksb(fjac,2,2,indx,p)
244 sig1=sig1+p(1)
245 sig2=sig2+p(2)
246 c1=log(-sig2*(1-sig1)/sig1/(sig2-1))/lupp
247 c2=log((sig2-1)/(1-sig1))
248 enddo
249 enddo
250! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
251! STEP 2.
252! Compute minimum surface pressure
253 ppred=ppre+dppre
254 spre=(ppre-pupp)/pdif
255 spred=(ppred-pupp)/pdif
256 rkpre=(log((spre-sig2)/(sig1-spre))-c2)/c1
257 rkpred=(log((spred-sig2)/(sig1-spred))-c2)/c1
258 pkpre=pdif*c1*(sig2-spre)*(spre-sig1)/(sig1-sig2)
259 pkkpre=pkpre*c1*(sig2+sig1-2*spre)/(sig1-sig2)
260 psigd=psig-dpsig
261 ssig=(psig-pupp)/pdif
262 ssigd=(psigd-pupp)/pdif
263 rksig=(log((ssig-sig2)/(sig1-ssig))-c2)/c1
264 rksigd=(log((ssigd-sig2)/(sig1-ssigd))-c2)/c1
265 pksig=pdif*c1*(sig2-ssig)*(ssig-sig1)/(sig1-sig2)
266 pkksig=pksig*c1*(sig2+sig1-2*ssig)/(sig1-sig2)
267 x2=0
268 a2=0
269 x2(1)=pkpre
270 a2(1,1)=1
271 a2(1,2)=rkpre
272 a2(1,3)=rkpre**2
273 x2(2)=pkkpre
274 a2(2,2)=1
275 a2(2,3)=2*rkpre
276 a2(3,1)=1
277 a2(3,2)=rkpred
278 a2(3,3)=rkpred**2
279 a2(3,4)=-1
280 a2(3,5)=-rkpred
281 a2(4,2)=1
282 a2(4,3)=2*rkpred
283 a2(4,5)=-1
284 a2(5,4)=-1
285 a2(5,5)=-rksigd
286 a2(5,6)=1
287 a2(5,7)=rksigd
288 a2(5,8)=rksigd**2
289 a2(6,5)=-1
290 a2(6,7)=1
291 a2(6,8)=2*rksigd
292 a2(7,6)=1
293 a2(7,7)=rksig
294 a2(7,8)=rksig**2
295 a2(7,9)=-pksig/pbot
296 a2(8,7)=1
297 a2(8,8)=2*rksig
298 a2(8,9)=-pkksig/pbot
299 x2(9)=ppre
300 a2(9,1)=(rkpre-rkpred)
301 a2(9,2)=(rkpre**2-rkpred**2)/2
302 a2(9,3)=(rkpre**3-rkpred**3)/3
303 a2(9,4)=(rkpred-rksigd)
304 a2(9,5)=(rkpred**2-rksigd**2)/2
305 a2(9,6)=(rksigd-rksig)
306 a2(9,7)=(rksigd**2-rksig**2)/2
307 a2(9,8)=(rksigd**3-rksig**3)/3
308 a2(9,9)=psig/pbot
309 call ludcmp(a2,9,9,indx2,d)
310 call lubksb(a2,9,9,indx2,x2)
311 pmin=x2(9)
312 p2sig=psig/pbot*pmin
313 p2sigd=p2sig &
314 +x2(6)*(rksigd-rksig) &
315 +x2(7)*(rksigd**2-rksig**2)/2 &
316 +x2(8)*(rksigd**3-rksig**3)/3
317 p2pred=p2sigd &
318 +x2(4)*(rkpred-rksigd) &
319 +x2(5)*(rkpred**2-rksigd**2)/2
320! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
321! STEP 3.
322 if(lupp.lt.levs) then
323 pkupp=pdif*c1*(sig2-0)*(0-sig1)/(sig1-sig2)
324 pkkupp=pkupp*c1*(sig2+sig1-2*0)/(sig1-sig2)
325 x3=0
326 a3=0
327 x3(1)=pkupp
328 a3(1,1)=pupp
329 x3(2)=pkkupp*pupp-pkupp**2
330 a3(2,2)=pupp**2
331 x3(3)=log((ptop+dptop)/pupp)
332 a3(3,1)=levs-1-lupp
333 a3(3,2)=(levs-1-lupp)**2/2
334 a3(3,3)=(levs-1-lupp)**3/3
335 call ludcmp(a3,3,3,indx3,d)
336 call lubksb(a3,3,3,indx3,x3)
337 endif
338! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
339! STEP 4.
340! Compute hybrid interface values
341 ak(0)=0
342 bk(0)=1
343 do k=1,levs-1
344 if(k.ge.lupp) then
345 p1=pupp*exp(x3(1)*(k-lupp)+x3(2)*(k-lupp)**2/2+x3(3)*(k-lupp)**3/3)
346 else
347 p1=(sig1*exp(c1*k+c2)+sig2)/(exp(c1*k+c2)+1)*pdif+pupp
348 endif
349 if(k.ge.rkpre) then
350 p2=p1
351 elseif(k.ge.rkpred) then
352 p2=p2pred+x2(1)*(k-rkpred) &
353 +x2(2)*(k**2-rkpred**2)/2 &
354 +x2(3)*(k**3-rkpred**3)/3
355 elseif(k.ge.rksigd) then
356 p2=p2sigd+x2(4)*(k-rksigd) &
357 +x2(5)*(k**2-rksigd**2)/2
358 elseif(k.ge.rksig) then
359 p2=p2sig+x2(6)*(k-rksig) &
360 +x2(7)*(k**2-rksig**2)/2 &
361 +x2(8)*(k**3-rksig**3)/3
362 else
363 p2=p1/pbot*pmin
364 endif
365 ak(k)=(p2*pbot-p1*pmin)/(pbot-pmin)
366 bk(k)=(p1-p2)/(pbot-pmin)
367 enddo
368 ak(levs)=ptop
369 bk(levs)=0
370end subroutine
subroutine lubksb(a, n, np, indx, b)
Lower and upper triangular back substitution.
subroutine ludcmp(a, n, np, indx, d)
This subprogram decomposes a matrix into a product of lower and upper triangular matrices.
subroutine vcoord_gen(levs, lupp, pbot, psig, ppre, pupp, ptop, dpbot, dpsig, dppre, dpupp, dptop, pmin, ak, bk)
This subprogram generates hybrid coordinate interface profiles from a few given parameters.