!-------------------------------------- LICENCE BEGIN ------------------------------------
!Environment Canada - Atmospheric Science and Technology License/Disclaimer,
! version 3; Last Modified: May 7, 2008.
!This is free but copyrighted software; you can use/redistribute/modify it under the terms
!of the Environment Canada - Atmospheric Science and Technology License/Disclaimer
!version 3 or (at your option) any later version that should be found at:
!http://collaboration.cmc.ec.gc.ca/science/rpn.comm/license.html
!
!This software is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY;
!without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
!See the above mentioned License/Disclaimer for more details.
!You should have received a copy of the License/Disclaimer along with this software;
!if not, you can write to: EC-RPN COMM Group, 2121 TransCanada, suite 500, Dorval (Quebec),
!CANADA, H9P 1J3; or send e-mail to service.rpn@ec.gc.ca
!-------------------------------------- LICENCE END --------------------------------------
***s/r e_intwind_offline - Wind components horizontal interpolation
*
#include "model_macros_f.h"
*
subroutine e_intwind_offline 2,7
implicit none
*
*author B. Bilodeau - spring 2008
*
*revision
* v3_31 - Bilodeau B. - initial version
* v3_31 - Lee V. - add 3DF pilot for Schm_offline_L
*
*object
* see above ID
*
*ARGUMENTS
*
*IMPLICITS
#include "e_option.cdk"
#include "e_fu.cdk"
#include "e_anal.cdk"
#include "e_grids.cdk"
#include "e_cdate.cdk"
#include "e_topo.cdk"
#include "dcst.cdk"
#include "grd.cdk"
#include "bmf.cdk"
#include "e_schm.cdk"
#include "pilot.cdk"
#include "e_grdc.cdk"
#include "hgc.cdk"
*
integer ezqkdef,ezdefset,ezsetopt,ezuvint,
$ fstinf,fstlir,fstprm,e_rdhint3
external ezqkdef,ezdefset,ezsetopt,ezuvint,
$ fstinf,fstlir,fstprm,e_rdhint3
*
integer i, j, k, src_gid, key1, key2, nic, njc, ni1, nj1,
$ nk1,nkc,err,iu,ju,iv,jv,nu,nv
integer nisu,nisv,njsu,njsv
integer nis,njs,niw,njw,iw,ie,jw,is,js,jn
integer ip2, ip3
integer dte, det, ipas, p1, p2, p3, g1, g2, g3, g4, bit,
$ dty, swa, lng, dlf, ubc, ex1, ex2, ex3
character*1 typ,grd
character*4 var,var_uu,var_vv
character*12 lab
logical must_interpo_L
real, dimension (:), allocatable :: uu,vv
real, dimension (:), allocatable :: w1,w2
real c1
real, dimension (:,:), allocatable :: uuf,vvf
real, dimension (:,:,:), allocatable :: uun,vvn
*
logical must_interpo_s
integer un_s,nic_s,njc_s,g1_s,g2_s,g3_s,kind,id_s
data un_s,nic_s,njc_s,g1_s,g2_s,g3_s,id_s /-1,-1,-1,-1,-1,-1,-1/
save un_s,nic_s,njc_s,g1_s,g2_s,g3_s,id_s,must_interpo_s
*
* ---------------------------------------------------------------------
*
write(6,1001)
if (Pil_bmf_L) then
nisu = niu
njsu = nju
nisv = niv
njsv = njv
nis = nifi
njs = njfi
else
nisu = E_Grdc_ni
njsu = E_Grdc_nj
nisv = E_Grdc_ni
njsv = E_Grdc_nj
nis = E_Grdc_ni
njs = E_Grdc_nj
allocate (uun(nisu,njsu,lv),vvn(nisv,njsv,lv))
endif
*
allocate (uu(nisu*njsu),vv(nisv*njsv))
* uuf and vvf are on phi grid
allocate (uuf(nis,njs),vvf(nis,njs))
*
nu = nisu*njsu
nv = nisv*njsv
*
ip2 = ip2a
ip3 = ip3a
*
key1=fstinf(e_fu_anal,nic,njc,nkc,datev,' ',na(1),ip2,ip3,
$ ' ','UU')
err= fstprm (key1, DTE, DET, IPAS, ni1, nj1, nk1, BIT, DTY, P1,
$ P2, P3, TYP, VAR, LAB, GRD, G1, G2, G3, G4, SWA,
$ LNG, DLF, UBC, EX1, EX2, EX3)
src_gid = ezqkdef (nic, njc, GRD, g1, g2, g3, g4, e_fu_anal)
*
* Should we interpolate the winds?
*
If (anal_hav(1).eq.0) then
must_interpo_L = .false.
print *,'NO interpolation required for winds'
if (nic.ne.nis.or.njc.ne.njs) then
print *,'ERROR: WINDS ARE NOT ON THE SAME GRID AS TT,HU'
print *,'nic=',nic,'.NE. nis=',nis,' njc=',njc,'.NE. njs=',njs
goto 55
endif
else
must_interpo_L = .true.
endif
*
allocate (w1(nic*njc),w2(nic*njc))
*
do 200 k=1,lv
*
ip2 = ip2a
ip3 = ip3a
key1 = fstlir (w1, e_fu_anal, iu, ju, nkc, datev, ' ',
$ na(k), ip2, ip3, ' ', 'UU')
key2 = fstlir (w2, e_fu_anal, iv, jv, nkc, datev, ' ',
$ na(k), ip2, ip3, ' ', 'VV')
if (key1.lt.0 .or. iu.ne.nic .or. ju.ne.njc ) then
write(6,*)'ERROR: UU NOT AVAILABLE,'
goto 55
endif
if (key1.lt.0 .or. iv.ne.nic .or. jv.ne.njc ) then
write(6,*)'ERROR: VV NOT AVAILABLE,'
goto 55
endif
*
* Horizontal Interpolation on PHI grid, if necessary
*
if (must_interpo_L) then
err = ezsetopt ('INTERP_DEGREE', 'CUBIC')
err = ezdefset ( dstf_gid, src_gid )
err = ezuvint ( uuf,vvf,w1,w2 )
write (6,*) '(S/R e_intwind_offline) Horizontal interpolation:UU,VV with ezuvint CUBIC '
else
*
write(6,'(A25,A30,A3,A12,3I5)')
$ "(S/R e_intwind_offline) ",
$ "NO horizontal interpolation on ","UU"," for ip123= ",na(k),ip2,ip3
write(6,'(A25,A30,A3,A12,3I5)')
$ "(S/R e_intwind_offline) ",
$ "NO horizontal interpolation on ","VV"," for ip123= ",na(k),ip2,ip3
do j=1,njs
do i=1,nis
uuf(i,j) = w1((j-1)*nis+i)
vvf(i,j) = w2((j-1)*nis+i)
end do
end do
endif
*
if (Pil_bmf_L) then
* copy uuf (phi grid) into uu
do j=1,njsu
do i=1,nisu
uu((j-1)*nisu+i) = uuf(i,j)
end do
end do
* copy vvf (phi grid) into vv
do j=1,njsv
do i=1,nisv
vv((j-1)*nisv+i) = vvf(i,j)
end do
end do
else
* copy and scale uuf (phi grid) into uu
do j=1,njsu
do i=1,nisu
uu((j-1)*nisu+i) = Dcst_knams_8*uuf(i,j)
end do
end do
* copy and scale vvf (phi grid) into vv
do j=1,njsv
do i=1,nisv
vv((j-1)*nisv+i) = Dcst_knams_8*vvf(i,j)
end do
end do
endif
*
if (Pil_bmf_L) then
call e_bmfsplitxy2
(uu,nisu,njsu,'UU ',k,lv,pniu,0,0,0)
call e_bmfsplitxy2
(vv,nisv,njsv,'VV ',k,lv,pni ,0,0,0)
else
call e_fill_3df
( uu,uun,nisu,njsu,lv,k,1.0,0.0)
call e_fill_3df
( vv,vvn,nisv,njsv,lv,k,1.0,0.0)
endif
*
200 continue
deallocate(uu,vv,uuf,vvf,w1,w2)
if (.not.Pil_bmf_L) then
call e_write_3df
( uun,nisu,njsu,lv,'UU ',unf_casc)
call e_write_3df
( vvn,nisv,njsv,lv,'VV ',unf_casc)
deallocate(uun,vvn)
return
endif
*
*
goto 99
55 call e_arret
( 'e_intwind' )
*
*
99 continue
call bmf_splitwrall ('AHAV',2,1,1,Bmf_time1,Bmf_time2,
$ 0,0,40,0,anal_hav)
call bmf_splitend
return
*
1001 format ('COMPUTE UU,VV (S/R E_INTWIND_OFFLINE)',/,25('+'))
120 format ('|',1x,a8,'|',1x,a2,' |',2(i7,' |'),i3,' |',1x,a3,
$ ' |',1x,a7,' |',1x,a16,'|')
* ---------------------------------------------------------------------
end