53 integer,
private,
allocatable :: interp_xi2z_idx (:,:,:,:)
54 real(RP),
private,
allocatable :: interp_xi2z_coef (:,:,:,:)
55 integer,
private,
allocatable :: interp_z2xi_idx (:,:,:,:)
56 real(RP),
private,
allocatable :: interp_z2xi_coef (:,:,:,:)
58 integer,
private,
allocatable :: interp_xi2p_idx (:,:,:,:)
59 real(RP),
private,
allocatable :: interp_xi2p_coef (:,:,:,:)
62 integer,
private,
allocatable :: interp_xih2zh_idx (:,:,:,:)
63 real(RP),
private,
allocatable :: interp_xih2zh_coef(:,:,:,:)
64 integer,
private,
allocatable :: interp_zh2xih_idx (:,:,:,:)
65 real(RP),
private,
allocatable :: interp_zh2xih_coef(:,:,:,:)
67 integer,
private,
allocatable :: interp_xih2p_idx (:,:,:,:)
68 real(RP),
private,
allocatable :: interp_xih2p_coef (:,:,:,:)
83 integer,
intent(in) :: KA, KS, KE
84 integer,
intent(in) :: IA, IS, IE
85 integer,
intent(in) :: JA, JS, JE
86 logical,
intent(in) :: TOPO_exist
87 real(RP),
intent(in) :: Xi ( ka)
88 real(RP),
intent(in) :: Xih(0:ka)
89 real(RP),
intent(in) :: Z ( ka,ia,ja)
90 real(RP),
intent(in) :: Zh (0:ka,ia,ja)
92 integer :: k, i, j, kk, kp
96 log_info(
"INTERP_VERT_setcoef",*)
'Setup' 97 log_info(
"INTERP_VERT_setcoef",*)
'No namelists.' 102 log_info(
"INTERP_VERT_setcoef",*)
'Topography exists & interpolation has meaning? : ',
interp_available 106 allocate( interp_xi2z_idx(ka,ia,ja,2) )
107 allocate( interp_xi2z_coef(ka,ia,ja,3) )
108 allocate( interp_z2xi_idx(ka,ia,ja,2) )
109 allocate( interp_z2xi_coef(ka,ia,ja,3) )
117 if ( xi(k) <= zh(ks-1,i,j) )
then 119 interp_xi2z_idx(k,i,j,1) = ks
120 interp_xi2z_idx(k,i,j,2) = ks
121 interp_xi2z_coef(k,i,j,1) = 0.0_rp
122 interp_xi2z_coef(k,i,j,2) = 0.0_rp
123 interp_xi2z_coef(k,i,j,3) = 1.0_rp
125 elseif( xi(k) <= z(ks,i,j) )
then 127 interp_xi2z_idx(k,i,j,1) = ks
128 interp_xi2z_idx(k,i,j,2) = ks
129 interp_xi2z_coef(k,i,j,1) = 0.0_rp
130 interp_xi2z_coef(k,i,j,2) = 1.0_rp
131 interp_xi2z_coef(k,i,j,3) = 0.0_rp
133 elseif( xi(k) > zh(ke,i,j) )
then 135 interp_xi2z_idx(k,i,j,1) = ke
136 interp_xi2z_idx(k,i,j,2) = ke
137 interp_xi2z_coef(k,i,j,1) = 0.0_rp
138 interp_xi2z_coef(k,i,j,2) = 0.0_rp
139 interp_xi2z_coef(k,i,j,3) = 1.0_rp
141 elseif( xi(k) > z(ke,i,j) )
then 143 interp_xi2z_idx(k,i,j,1) = ke
144 interp_xi2z_idx(k,i,j,2) = ke
145 interp_xi2z_coef(k,i,j,1) = 1.0_rp
146 interp_xi2z_coef(k,i,j,2) = 0.0_rp
147 interp_xi2z_coef(k,i,j,3) = 0.0_rp
153 if( xi(k) <= z(kk,i,j) )
exit 156 interp_xi2z_idx(k,i,j,1) = kp - 1
157 interp_xi2z_idx(k,i,j,2) = kp
158 interp_xi2z_coef(k,i,j,1) = ( z(kp,i,j) - xi(k) ) &
159 / ( z(kp,i,j) - z(kp-1,i,j) )
160 interp_xi2z_coef(k,i,j,2) = ( xi(k) - z(kp-1,i,j) ) &
161 / ( z(kp,i,j) - z(kp-1,i,j) )
162 interp_xi2z_coef(k,i,j,3) = 0.0_rp
175 if ( z(k,i,j) < xih(ks-1) )
then 177 interp_z2xi_idx(k,i,j,1) = ks
178 interp_z2xi_idx(k,i,j,2) = ks
179 interp_z2xi_coef(k,i,j,1) = 0.0_rp
180 interp_z2xi_coef(k,i,j,2) = 1.0_rp
183 elseif( z(k,i,j) <= xi(ks) )
then 185 interp_z2xi_idx(k,i,j,1) = ks
186 interp_z2xi_idx(k,i,j,2) = ks
187 interp_z2xi_coef(k,i,j,1) = 0.0_rp
188 interp_z2xi_coef(k,i,j,2) = 1.0_rp
189 interp_z2xi_coef(k,i,j,3) = 0.0_rp
191 elseif( z(k,i,j) > xih(ke) )
then 193 interp_z2xi_idx(k,i,j,1) = ke
194 interp_z2xi_idx(k,i,j,2) = ke
195 interp_z2xi_coef(k,i,j,1) = 0.0_rp
196 interp_z2xi_coef(k,i,j,2) = 0.0_rp
197 interp_z2xi_coef(k,i,j,3) = 1.0_rp
199 elseif( z(k,i,j) > xi(ke) )
then 201 interp_z2xi_idx(k,i,j,1) = ke
202 interp_z2xi_idx(k,i,j,2) = ke
203 interp_z2xi_coef(k,i,j,1) = 1.0_rp
204 interp_z2xi_coef(k,i,j,2) = 0.0_rp
205 interp_z2xi_coef(k,i,j,3) = 0.0_rp
211 if( z(k,i,j) <= xi(kk) )
exit 214 interp_z2xi_idx(k,i,j,1) = kp - 1
215 interp_z2xi_idx(k,i,j,2) = kp
216 interp_z2xi_coef(k,i,j,1) = ( xi(kp) - z(k,i,j) ) &
217 / ( xi(kp) - xi(kp-1) )
218 interp_z2xi_coef(k,i,j,2) = ( z(k,i,j) - xi(kp-1) ) &
219 / ( xi(kp) - xi(kp-1) )
220 interp_z2xi_coef(k,i,j,3) = 0.0_rp
232 interp_xi2z_idx( 1:ks-1,i,j,1) = ks
233 interp_xi2z_idx( 1:ks-1,i,j,2) = ks
234 interp_xi2z_coef( 1:ks-1,i,j,1) = 0.0_rp
235 interp_xi2z_coef( 1:ks-1,i,j,2) = 0.0_rp
236 interp_xi2z_coef( 1:ks-1,i,j,3) = 1.0_rp
238 interp_xi2z_idx(ke+1:ka,i,j,1) = ke
239 interp_xi2z_idx(ke+1:ka,i,j,2) = ke
240 interp_xi2z_coef(ke+1:ka,i,j,1) = 0.0_rp
241 interp_xi2z_coef(ke+1:ka,i,j,2) = 0.0_rp
242 interp_xi2z_coef(ke+1:ka,i,j,3) = 1.0_rp
244 interp_z2xi_idx( 1:ks-1,i,j,1) = ks
245 interp_z2xi_idx( 1:ks-1,i,j,2) = ks
246 interp_z2xi_coef( 1:ks-1,i,j,1) = 0.0_rp
247 interp_z2xi_coef( 1:ks-1,i,j,2) = 0.0_rp
248 interp_z2xi_coef( 1:ks-1,i,j,3) = 1.0_rp
250 interp_z2xi_idx(ke+1:ka,i,j,1) = ke
251 interp_z2xi_idx(ke+1:ka,i,j,2) = ke
252 interp_z2xi_coef(ke+1:ka,i,j,1) = 0.0_rp
253 interp_z2xi_coef(ke+1:ka,i,j,2) = 0.0_rp
254 interp_z2xi_coef(ke+1:ka,i,j,3) = 1.0_rp
261 allocate( interp_xih2zh_idx(0:ka,ia,ja,2) )
262 allocate( interp_xih2zh_coef(0:ka,ia,ja,3) )
263 allocate( interp_zh2xih_idx(0:ka,ia,ja,2) )
264 allocate( interp_zh2xih_coef(0:ka,ia,ja,3) )
272 if ( xih(k) < zh(ks-1,i,j) )
then 274 interp_xih2zh_idx(k,i,j,1) = ks
275 interp_xih2zh_idx(k,i,j,2) = ks
276 interp_xih2zh_coef(k,i,j,1) = 0.0_rp
277 interp_xih2zh_coef(k,i,j,2) = 0.0_rp
278 interp_xih2zh_coef(k,i,j,3) = 1.0_rp
280 elseif( xih(k) > zh(ke,i,j) )
then 282 interp_xih2zh_idx(k,i,j,1) = ke
283 interp_xih2zh_idx(k,i,j,2) = ke
284 interp_xih2zh_coef(k,i,j,1) = 0.0_rp
285 interp_xih2zh_coef(k,i,j,2) = 0.0_rp
286 interp_xih2zh_coef(k,i,j,3) = 1.0_rp
292 if( xih(k) <= zh(kk,i,j) )
exit 295 interp_xih2zh_idx(k,i,j,1) = kp - 1
296 interp_xih2zh_idx(k,i,j,2) = kp
297 interp_xih2zh_coef(k,i,j,1) = ( zh(kp,i,j) - xih(k) ) &
298 / ( zh(kp,i,j) - zh(kp-1,i,j) )
299 interp_xih2zh_coef(k,i,j,2) = ( xih(k) - zh(kp-1,i,j) ) &
300 / ( zh(kp,i,j) - zh(kp-1,i,j) )
301 interp_xih2zh_coef(k,i,j,3) = 0.0_rp
314 if ( zh(k,i,j) <= xih(ks-1) )
then 316 interp_zh2xih_idx(k,i,j,1) = ks
317 interp_zh2xih_idx(k,i,j,2) = ks
318 interp_zh2xih_coef(k,i,j,1) = 0.0_rp
319 interp_zh2xih_coef(k,i,j,2) = 0.0_rp
320 interp_zh2xih_coef(k,i,j,3) = 1.0_rp
322 elseif( zh(k,i,j) > xih(ke) )
then 324 interp_zh2xih_idx(k,i,j,1) = ke
325 interp_zh2xih_idx(k,i,j,2) = ke
326 interp_zh2xih_coef(k,i,j,1) = 0.0_rp
327 interp_zh2xih_coef(k,i,j,2) = 0.0_rp
328 interp_zh2xih_coef(k,i,j,3) = 1.0_rp
334 if( zh(k,i,j) <= xih(kk) )
exit 337 interp_zh2xih_idx(k,i,j,1) = kp - 1
338 interp_zh2xih_idx(k,i,j,2) = kp
339 interp_zh2xih_coef(k,i,j,1) = ( xih(kp) - zh(k,i,j) ) &
340 / ( xih(kp) - xih(kp-1) )
341 interp_zh2xih_coef(k,i,j,2) = ( zh(k,i,j) - xih(kp-1) ) &
342 / ( xih(kp) - xih(kp-1) )
343 interp_zh2xih_coef(k,i,j,3) = 0.0_rp
355 interp_xih2zh_idx( 1:ks-2,i,j,1) = ks
356 interp_xih2zh_idx( 1:ks-2,i,j,2) = ks
357 interp_xih2zh_coef( 1:ks-2,i,j,1) = 0.0_rp
358 interp_xih2zh_coef( 1:ks-2,i,j,2) = 0.0_rp
359 interp_xih2zh_coef( 1:ks-2,i,j,3) = 1.0_rp
361 interp_xih2zh_idx(ke+1:ka,i,j,1) = ke
362 interp_xih2zh_idx(ke+1:ka,i,j,2) = ke
363 interp_xih2zh_coef(ke+1:ka,i,j,1) = 0.0_rp
364 interp_xih2zh_coef(ke+1:ka,i,j,2) = 0.0_rp
365 interp_xih2zh_coef(ke+1:ka,i,j,3) = 1.0_rp
367 interp_zh2xih_idx( 1:ks-2,i,j,1) = ks
368 interp_zh2xih_idx( 1:ks-2,i,j,2) = ks
369 interp_zh2xih_coef( 1:ks-2,i,j,1) = 0.0_rp
370 interp_zh2xih_coef( 1:ks-2,i,j,2) = 0.0_rp
371 interp_zh2xih_coef( 1:ks-2,i,j,3) = 1.0_rp
373 interp_zh2xih_idx(ke+1:ka,i,j,1) = ke
374 interp_zh2xih_idx(ke+1:ka,i,j,2) = ke
375 interp_zh2xih_coef(ke+1:ka,i,j,1) = 0.0_rp
376 interp_zh2xih_coef(ke+1:ka,i,j,2) = 0.0_rp
377 interp_zh2xih_coef(ke+1:ka,i,j,3) = 1.0_rp
396 matrix_solver_tridiagonal
399 integer,
intent(in) :: KA, KS, KE
400 integer,
intent(in) :: IA, IS, IE
401 integer,
intent(in) :: JA, JS, JE
402 real(RP),
intent(in) :: Xi (ka)
403 real(RP),
intent(in) :: Z (ka,ia,ja)
404 real(RP),
intent(in) :: var (ka,ia,ja)
405 real(RP),
intent(out) :: var_Z(ka,ia,ja)
411 real(RP) :: c1, c2, c3, d
414 integer :: k, i, j, kk
419 if ( kmax == 2 )
then 429 var_z(k,i,j) = interp_xi2z_coef(k,i,j,1) * var(interp_xi2z_idx(k,i,j,1),i,j) &
430 + interp_xi2z_coef(k,i,j,2) * var(interp_xi2z_idx(k,i,j,2),i,j) &
448 fdz(k) = z(k+1,i,j) - z(k,i,j)
451 md(ks+1) = 2.0 * ( fdz(ks) + fdz(ks+1) ) + fdz(ks)
453 md(k) = 2.0 * ( fdz(k-1) + fdz(k) )
455 md(ke-1) = 2.0 * ( fdz(ke-2) + fdz(ke-1) ) + fdz(ke-1)
458 v(k) = ( var(k+1,i,j) - var(k ,i,j) ) / fdz(k ) &
459 - ( var(k ,i,j) - var(k-1,i,j) ) / fdz(k-1)
462 call matrix_solver_tridiagonal( ka, ks+1, ke-1, &
473 kk = min(interp_xi2z_idx(k,i,j,1),ke-1)
475 c3 = ( u(kk+1) - u(kk) ) / fdz(kk)
477 c1 = ( var(kk+1,i,j) - var(kk,i,j) ) / fdz(kk) - ( u(kk) * 2.0_rp + u(kk+1) ) * fdz(kk)
478 d = xi(k) - z(kk,i,j)
480 var_z(k,i,j) = ( interp_xi2z_coef(k,i,j,3) ) *
const_undef &
481 + ( 1.0_rp-interp_xi2z_coef(k,i,j,3) ) * ( ( ( c3 * d + c2 ) * d + c1 ) * d + var(kk,i,j) )
504 matrix_solver_tridiagonal
507 integer,
intent(in) :: KA, KS, KE
508 integer,
intent(in) :: IA, IS, IE
509 integer,
intent(in) :: JA, JS, JE
510 real(RP),
intent(in) :: Z (ka,ia,ja)
511 real(RP),
intent(in) :: Xi (ka)
512 real(RP),
intent(in) :: var (ka,ia,ja)
513 real(RP),
intent(out) :: var_Xi(ka,ia,ja)
519 real(RP) :: c1, c2, c3, d
522 integer :: k, i, j, kk
527 if ( kmax <= 2 )
then 537 var_xi(k,i,j) = interp_z2xi_coef(k,i,j,1) * var(interp_z2xi_idx(k,i,j,1),i,j) &
538 + interp_z2xi_coef(k,i,j,2) * var(interp_z2xi_idx(k,i,j,2),i,j) &
547 fdz(k) = xi(k+1) - xi(k)
550 md(ks+1) = 2.0 * ( fdz(ks) + fdz(ks+1) ) + fdz(ks)
552 md(k) = 2.0 * ( fdz(k-1) + fdz(k) )
554 md(ke-1) = 2.0 * ( fdz(ke-2) + fdz(ke-1) ) + fdz(ke-1)
566 v(k) = ( var(k+1,i,j) - var(k ,i,j) ) / fdz(k ) &
567 - ( var(k ,i,j) - var(k-1,i,j) ) / fdz(k-1)
570 call matrix_solver_tridiagonal( ka, ks+1, ke-1, &
581 kk = min(interp_z2xi_idx(k,i,j,1),ke-1)
583 c3 = ( u(kk+1) - u(kk) ) / fdz(kk)
585 c1 = ( var(kk+1,i,j) - var(kk,i,j) ) / fdz(kk) - ( u(kk) * 2.0_rp + u(kk+1) ) * fdz(kk)
586 d = z(k,i,j) - xi(kk)
588 var_xi(k,i,j) = ( interp_z2xi_coef(k,i,j,3) ) *
const_undef &
589 + ( 1.0_rp-interp_z2xi_coef(k,i,j,3) ) * ( ( ( c3 * d + c2 ) * d + c1 ) * d + var(kk,i,j) )
612 matrix_solver_tridiagonal
615 integer,
intent(in) :: KA, KS, KE
616 integer,
intent(in) :: IA, IS, IE
617 integer,
intent(in) :: JA, JS, JE
618 real(RP),
intent(in) :: Xih (0:ka)
619 real(RP),
intent(in) :: Zh (0:ka,ia,ja)
620 real(RP),
intent(in) :: var (ka,ia,ja)
621 real(RP),
intent(out) :: var_Z(ka,ia,ja)
627 real(RP) :: c1, c2, c3, d
630 integer :: k, i, j, kk
635 if ( kmax == 2 )
then 645 var_z(k,i,j) = interp_xih2zh_coef(k,i,j,1) * var(interp_xih2zh_idx(k,i,j,1),i,j) &
646 + interp_xih2zh_coef(k,i,j,2) * var(interp_xih2zh_idx(k,i,j,2),i,j) &
664 cdz(k) = zh(k,i,j) - zh(k-1,i,j)
667 md(ks) = 2.0 * ( cdz(ks) + cdz(ks+1) ) + cdz(ks)
669 md(k) = 2.0 * ( cdz(k) + cdz(k+1) )
671 md(ke-1) = 2.0 * ( cdz(ke-1) + cdz(ke) ) + cdz(ke)
674 v(k) = ( var(k+1,i,j) - var(k ,i,j) ) / cdz(k+1) &
675 - ( var(k ,i,j) - var(k-1,i,j) ) / cdz(k )
678 call matrix_solver_tridiagonal( kmax-1, 1, kmax-1, &
689 kk = min(interp_xih2zh_idx(k,i,j,1),ke-1)
691 c3 = ( u(kk+1) - u(kk) ) / cdz(kk+1)
693 c1 = ( var(kk+1,i,j) - var(kk,i,j) ) / cdz(kk+1) - ( u(kk) * 2.0_rp + u(kk+1) ) * cdz(kk+1)
694 d = xih(k) - zh(kk,i,j)
696 var_z(k,i,j) = ( interp_xih2zh_coef(k,i,j,3) ) *
const_undef &
697 + ( 1.0_rp-interp_xih2zh_coef(k,i,j,3) ) * ( ( ( c3 * d + c2 ) * d + c1 ) * d + var(kk,i,j) )
720 matrix_solver_tridiagonal
723 integer,
intent(in) :: KA, KS, KE
724 integer,
intent(in) :: IA, IS, IE
725 integer,
intent(in) :: JA, JS, JE
726 real(RP),
intent(in) :: Zh (0:ka,ia,ja)
727 real(RP),
intent(in) :: Xih (0:ka)
728 real(RP),
intent(in) :: var (ka,ia,ja)
729 real(RP),
intent(out) :: var_Xi(ka,ia,ja)
735 real(RP) :: c1, c2, c3, d
738 integer :: k, i, j, kk
743 if ( kmax == 2 )
then 753 var_xi(k,i,j) = interp_zh2xih_coef(k,i,j,1) * var(interp_zh2xih_idx(k,i,j,1),i,j) &
754 + interp_zh2xih_coef(k,i,j,2) * var(interp_zh2xih_idx(k,i,j,2),i,j) &
763 cdz(k) = xih(k) - xih(k-1)
766 md(ks) = 2.0 * ( cdz(ks) + cdz(ks+1) ) + cdz(ks)
768 md(k) = 2.0 * ( cdz(k) + cdz(k+1) )
770 md(ke-1) = 2.0 * ( cdz(ke-1) + cdz(ke) ) + cdz(ke)
782 v(k) = ( var(k+1,i,j) - var(k ,i,j) ) / cdz(k+1) &
783 - ( var(k ,i,j) - var(k-1,i,j) ) / cdz(k )
786 call matrix_solver_tridiagonal( kmax-1, 1, kmax-1, &
797 kk = min(interp_zh2xih_idx(k,i,j,1),ke-1)
799 c3 = ( u(kk+1) - u(kk) ) / cdz(kk+1)
801 c1 = ( var(kk+1,i,j) - var(kk,i,j) ) / cdz(kk+1) - ( u(kk) * 2.0_rp + u(kk+1) ) * cdz(kk+1)
802 d = zh(k,i,j) - xih(kk)
803 var_xi(k,i,j) = interp_zh2xih_coef(k,i,j,3) *
const_undef &
804 + ( 1.0_rp - interp_zh2xih_coef(k,i,j,3) ) &
805 * ( ( ( c3 * d + c2 ) * d + c1 ) * d + var(kk,i,j) )
822 integer,
intent(in) :: Kpres
823 integer,
intent(in) :: IA
824 integer,
intent(in) :: JA
827 allocate( interp_xi2p_idx(kpres,ia,ja,2) )
828 allocate( interp_xi2p_coef(kpres,ia,ja,3) )
829 allocate( interp_xih2p_idx(kpres,ia,ja,2) )
830 allocate( interp_xih2p_coef(kpres,ia,ja,3) )
847 integer,
intent(in) :: Kpres
848 integer,
intent(in) :: KA, KS, KE
849 integer,
intent(in) :: IA, IS, IE
850 integer,
intent(in) :: JA, JS, JE
851 real(RP),
intent(in) :: PRES ( ka,ia,ja)
852 real(RP),
intent(in) :: PRESh (0:ka,ia,ja)
853 real(RP),
intent(in) :: SFC_PRES( ia,ja)
854 real(RP),
intent(in) :: Paxis (kpres)
856 real(RP) :: LnPRES ( ka,ia,ja)
857 real(RP) :: LnPRESh (0:ka,ia,ja)
858 real(RP) :: LnSFC_PRES( ia,ja)
859 real(RP) :: LnPaxis (kpres)
861 integer :: k, i, j, kk, kp
867 lnpres(k,i,j) = log( pres(k,i,j) )
875 lnpresh(k,i,j) = log( presh(k,i,j) )
882 lnsfc_pres(i,j) = log( sfc_pres(i,j) )
886 lnpaxis(:) = log( paxis(:) )
893 if ( lnpaxis(k) >= lnsfc_pres(i,j) )
then 895 interp_xi2p_idx(k,i,j,1) = ks
896 interp_xi2p_idx(k,i,j,2) = ks
897 interp_xi2p_coef(k,i,j,1) = 0.0_rp
898 interp_xi2p_coef(k,i,j,2) = 0.0_rp
899 interp_xi2p_coef(k,i,j,3) = 1.0_rp
901 elseif( lnpaxis(k) >= lnpres(ks,i,j) )
then 903 interp_xi2p_idx(k,i,j,1) = ks
904 interp_xi2p_idx(k,i,j,2) = ks
905 interp_xi2p_coef(k,i,j,1) = 0.0_rp
906 interp_xi2p_coef(k,i,j,2) = 1.0_rp
907 interp_xi2p_coef(k,i,j,3) = 0.0_rp
909 elseif( lnpaxis(k) < lnpres(ke,i,j) )
then 911 interp_xi2p_idx(k,i,j,1) = ke
912 interp_xi2p_idx(k,i,j,2) = ke
913 interp_xi2p_coef(k,i,j,1) = 0.0_rp
914 interp_xi2p_coef(k,i,j,2) = 0.0_rp
915 interp_xi2p_coef(k,i,j,3) = 1.0_rp
921 if( lnpaxis(k) >= lnpres(kk,i,j) )
exit 924 interp_xi2p_idx(k,i,j,1) = kp - 1
925 interp_xi2p_idx(k,i,j,2) = kp
926 interp_xi2p_coef(k,i,j,1) = ( lnpaxis(k) - lnpres(kp,i,j) ) &
927 / ( lnpres(kp-1,i,j) - lnpres(kp,i,j) )
928 interp_xi2p_coef(k,i,j,2) = ( lnpres(kp-1,i,j) - lnpaxis(k) ) &
929 / ( lnpres(kp-1,i,j) - lnpres(kp,i,j) )
930 interp_xi2p_coef(k,i,j,3) = 0.0_rp
942 if ( lnpaxis(k) > lnsfc_pres(i,j) )
then 944 interp_xih2p_idx(k,i,j,1) = ks
945 interp_xih2p_idx(k,i,j,2) = ks
946 interp_xih2p_coef(k,i,j,1) = 0.0_rp
947 interp_xih2p_coef(k,i,j,2) = 0.0_rp
948 interp_xih2p_coef(k,i,j,3) = 1.0_rp
950 elseif( lnpaxis(k) < lnpresh(ke,i,j) )
then 952 interp_xih2p_idx(k,i,j,1) = ke
953 interp_xih2p_idx(k,i,j,2) = ke
954 interp_xih2p_coef(k,i,j,1) = 0.0_rp
955 interp_xih2p_coef(k,i,j,2) = 0.0_rp
956 interp_xih2p_coef(k,i,j,3) = 1.0_rp
962 if( lnpaxis(k) >= lnpresh(kk,i,j) )
exit 965 interp_xih2p_idx(k,i,j,1) = kp - 1
966 interp_xih2p_idx(k,i,j,2) = kp
967 interp_xih2p_coef(k,i,j,1) = ( lnpaxis(k) - lnpresh(kp,i,j) ) &
968 / ( lnpresh(kp-1,i,j) - lnpresh(kp,i,j) )
969 interp_xih2p_coef(k,i,j,2) = ( lnpresh(kp-1,i,j) - lnpaxis(k) ) &
970 / ( lnpresh(kp-1,i,j) - lnpresh(kp,i,j) )
971 interp_xih2p_coef(k,i,j,3) = 0.0_rp
993 integer,
intent(in) :: Kpres
994 integer,
intent(in) :: KA
995 integer,
intent(in) :: IA, IS, IE
996 integer,
intent(in) :: JA, JS, JE
997 real(RP),
intent(in) :: var (ka ,ia,ja)
998 real(RP),
intent(out) :: var_P(kpres,ia,ja)
1011 var_p(k,i,j) = interp_xi2p_coef(k,i,j,1) * var(interp_xi2p_idx(k,i,j,1),i,j) &
1012 + interp_xi2p_coef(k,i,j,2) * var(interp_xi2p_idx(k,i,j,2),i,j) &
1033 integer,
intent(in) :: Kpres
1034 integer,
intent(in) :: KA
1035 integer,
intent(in) :: IA, IS, IE
1036 integer,
intent(in) :: JA, JS, JE
1037 real(RP),
intent(in) :: var (ka ,ia,ja)
1038 real(RP),
intent(out) :: var_P(kpres,ia,ja)
1051 var_p(k,i,j) = interp_xih2p_coef(k,i,j,1) * var(interp_xih2p_idx(k,i,j,1),i,j) &
1052 + interp_xih2p_coef(k,i,j,2) * var(interp_xih2p_idx(k,i,j,2),i,j) &
subroutine, public interp_vert_xi2p(Kpres, KA, IA, IS, IE, JA, JS, JE, var, var_P)
subroutine, public interp_vert_z2xi(KA, KS, KE, IA, IS, IE, JA, JS, JE, Z, Xi, var, var_Xi)
logical, public interp_available
topography exists & vertical interpolation has meaning?
module INTERPOLATION vertical
subroutine, public interp_vert_xi2z(KA, KS, KE, IA, IS, IE, JA, JS, JE, Xi, Z, var, var_Z)
real(rp), public const_undef
module atmosphere / grid / cartesC index
subroutine, public interp_vert_alloc_pres(Kpres, IA, JA)
Setup.
subroutine, public interp_vert_xih2zh(KA, KS, KE, IA, IS, IE, JA, JS, JE, Xih, Zh, var, var_Z)
subroutine, public interp_vert_setcoef_pres(Kpres, KA, KS, KE, IA, IS, IE, JA, JS, JE, PRES, PRESh, SFC_PRES, Paxis)
subroutine, public interp_vert_zh2xih(KA, KS, KE, IA, IS, IE, JA, JS, JE, Zh, Xih, var, var_Xi)
subroutine, public interp_vert_xih2p(Kpres, KA, IA, IS, IE, JA, JS, JE, var, var_P)
subroutine, public interp_vert_setcoef(KA, KS, KE, IA, IS, IE, JA, JS, JE, TOPO_exist, Xi, Xih, Z, Zh)
Setup.