!-------------------------------------- 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 --------------------------------------
*** casc_hvi - take incoming lo-res data, interpolate to hi-res data
* and place into the pieces of the BCS target grid as
* described in bcs_did
*
#include "model_macros_f.h"
*
subroutine casc_hvi (trname_a, 4,35
$ xpqd,ypqd,xpud,ypvd,xpqs,ypqs,xpus,ypvs,
$ uu1,uu2,vv1,vv2,tt1,tt2,psd1,psd2,
$ pip1,pip2,fip1,fip2,td1,td2,fi1,fi2,
$ qq1,qq2,tp1,tp2,tr1,ww1,ww2,mu1,mu2,ss1,ss2,ib1,ib2,
$ topo1,topo2,topu1,topu2,topv1,topv2,
$ uun,vvn,psdn,ttn,tpn,tdn,fin,qqn,ssn,fipn,pipn,wwn,mun,trn,
$ n1,n2,n3,n4,d1,d2,lnk,nid,njd,
$ nis,njs,nka,nvar,ntra,b1,b2)
implicit none
*
character*8 trname_a(ntra)
logical b1,b2
integer n1,n2,n3,n4,d1,d2,lnk,nid,njd,nis,njs,nka
integer ii,jj,jjj,kk,nvar,ntra,ib1,ib2,dd1,dd2
real*8 xpqd(*),ypqd(*),xpud(*),ypvd(*)
real*8 xpqs(*),ypqs(*),xpus(*),ypvs(*)
real uun(*),vvn(*),psdn(*),ttn(*),tpn(*),tdn(*),fin(nis*njs,nka),
$ qqn(*),ssn(*),fipn(*),
$ pipn(nis*njs,nka),wwn(*),mun(*),trn(*),
$ uu1(n1:n2,n3:n4,*),uu2(n1:n2,n3:n4,*),
$ tt1(n1:n2,n3:n4,*),tt2(n1:n2,n3:n4,*),
$ vv1(n1:n2,n3:n4,*),vv2(n1:n2,n3:n4,*),
$ psd1(n1:n2,n3:n4,*),psd2(n1:n2,n3:n4,*),
$ pip1(n1:n2,n3:n4,*),pip2(n1:n2,n3:n4,*),
$ fip1(n1:n2,n3:n4,*),fip2(n1:n2,n3:n4,*),
$ td1(n1:n2,n3:n4,*),td2(n1:n2,n3:n4,*),
$ fi1(n1:n2,n3:n4,*),fi2(n1:n2,n3:n4,*),
$ qq1(n1:n2,n3:n4,*),qq2(n1:n2,n3:n4,*),
$ tp1(n1:n2,n3:n4,*),tp2(n1:n2,n3:n4,*),
$ tr1(*),
$ ww1(n1:n2,n3:n4,*),ww2(n1:n2,n3:n4,*),
$ mu1(n1:n2,n3:n4,*),mu2(n1:n2,n3:n4,*),
$ ss1(n1:n2,n3:n4,*),ss2(n1:n2,n3:n4,*),
$ topo1(n1:n2,n3:n4),topo2(n1:n2,n3:n4),
$ topu1(n1:n2,n3:n4),topu2(n1:n2,n3:n4),
$ topv1(n1:n2,n3:n4),topv2(n1:n2,n3:n4)
*
*author M. Desgagne 2001 (MC2)
*
*revision
* v3_30 - Lee V. - initial version for GEMLAM
* v3_31 - Lee V. - correction in interpolation of firu,firv
#include "glb_ld.cdk"
#include "bcsdim.cdk"
#include "bcsgrds.cdk"
#include "bcsmem.cdk"
#include "lam.cdk"
#include "geomg.cdk"
#include "schm.cdk"
#include "pres.cdk"
#include "tr3d.cdk"
#include "itf_phy_busind.cdk"
*
integer i,j,k,n,ngd,nga,err,cnt,id
logical Vertint_L
integer, dimension (:), allocatable :: idx,idu,idy
real , dimension (:), allocatable :: wlnph,ana_p0,ana_p0u,ana_p0v
real, dimension (:,:), allocatable :: uur,vvr,psdr,ttr,tpr,tdr,
$ fir, qqr,ssr,fipr,pipr,wwr,mur
real, dimension (:,:,:), allocatable :: trr,uu3,vv3,psd3,tt3,tp3,
$ td3,fi3,qq3,ss3,fip3,pip3,ww3,mu3,tr3
real, dimension (:,:), allocatable ::
$ gz_temp,ttru,ttrv,firu,firv,ps,psu,psv,
$ topo_temp,topu_temp,topv_temp
real*8, dimension (: ), allocatable ::
$ cxa,cxb,cxc,cxd,cua,cub,cuc,cud,cya,cyb,cyc,cyd
*
*-----------------------------------------------------------------------
*
ngd = nid * njd
nga = nis * njs
if (ngd.le.0) return
*
allocate ( idx(nid), idu(max(nid,njd)), idy(njd) )
allocate ( cxa(nid),cxb(nid),cxc(nid),cxd(nid),
$ cua(max(nid,njd)),cub(max(nid,njd)),
$ cuc(max(nid,njd)),cud(max(nid,njd)),
$ cya(njd),cyb(njd),cyc(njd),cyd(njd))
*
call grid_to_grid_coef(xpqd,nid,xpqs,nis,idx,cxa,cxb,cxc,cxd,
$ Lam_hint_S)
call grid_to_grid_coef(ypqd,njd,ypqs,njs,idy,cya,cyb,cyc,cyd,
$ Lam_hint_S)
*
allocate (uur(ngd,nka),vvr(ngd,nka),psdr(ngd,nka),
$ ttr(ngd,nka),tpr(ngd,nka),tdr(ngd,nka), fir(ngd,nka),
$ qqr(ngd,nka),ssr(ngd,1),fipr(ngd,nka),pipr(ngd,nka),
$ wwr(ngd,nka),mur(ngd,nka),wlnph(nga),
$ ana_p0(ngd), ana_p0u(ngd),ana_p0v(ngd),
$ firu(ngd,nka),firv(ngd,nka),ttru(ngd,nka),ttrv(ngd,nka),
$ trr(ngd,nka,ntra))
*
* Perform horizontal interpolations
*
if (nvar.eq.5)then
Vertint_L=.true.
call hinterpo
( ttr,nid,njd, ttn,nis,njs,nka,
$ idx,idy,cxa,cxb,cxc,cxd,cya,cyb,cyc,cyd,Lam_hint_S)
call hinterpo
( fir(1,nka),nid,njd, fin( 1,nka),nis,njs,1,
$ idx,idy,cxa,cxb,cxc,cxd,cya,cyb,cyc,cyd,Lam_hint_S)
else
Vertint_L=.false.
call hinterpo
(ttr,nid,njd,ttn,nis,njs,nka,
$ idx,idy,cxa,cxb,cxc,cxd,cya,cyb,cyc,cyd,Lam_hint_S)
call hinterpo
(fir,nid,njd,fin,nis,njs,nka,
$ idx,idy,cxa,cxb,cxc,cxd,cya,cyb,cyc,cyd,Lam_hint_S)
call hinterpo
(qqr,nid,njd,qqn,nis,njs,nka,
$ idx,idy,cxa,cxb,cxc,cxd,cya,cyb,cyc,cyd,Lam_hint_S)
call hinterpo
(pipr,nid,njd,pipn,nis,njs,nka,
$ idx,idy,cxa,cxb,cxc,cxd,cya,cyb,cyc,cyd,Lam_hint_S)
call hinterpo
(tpr,nid,njd,tpn,nis,njs,nka,
$ idx,idy,cxa,cxb,cxc,cxd,cya,cyb,cyc,cyd,Lam_hint_S)
call hinterpo
(ssr,nid,njd,ssn,nis,njs,1,
$ idx,idy,cxa,cxb,cxc,cxd,cya,cyb,cyc,cyd,Lam_hint_S)
call hinterpo
(fipr,nid,njd,fipn,nis,njs,nka,
$ idx,idy,cxa,cxb,cxc,cxd,cya,cyb,cyc,cyd,Lam_hint_S)
call hinterpo
(psdr,nid,njd,psdn,nis,njs,nka,
$ idx,idy,cxa,cxb,cxc,cxd,cya,cyb,cyc,cyd,Lam_hint_S)
call hinterpo
(tdr,nid,njd,tdn,nis,njs,nka,
$ idx,idy,cxa,cxb,cxc,cxd,cya,cyb,cyc,cyd,Lam_hint_S)
If (.not.Schm_hydro_L.and.nvar.gt.11) then
call hinterpo
(wwr,nid,njd,wwn,nis,njs,nka,
$ idx,idy,cxa,cxb,cxc,cxd,cya,cyb,cyc,cyd,Lam_hint_S)
call hinterpo
(mur,nid,njd,mun,nis,njs,nka,
$ idx,idy,cxa,cxb,cxc,cxd,cya,cyb,cyc,cyd,Lam_hint_S)
endif
endif
* Compute wlnph, surface pressure from analysis
do i=1,nis*njs
wlnph(i) = ana_z(nka)+pipn(i,nka)
enddo
call hinterpo
(ana_p0,nid,njd,wlnph,nis,njs,1,
$ idx,idy,cxa,cxb,cxc,cxd,cya,cyb,cyc,cyd,Lam_hint_S)
*
do k=1,ntra
if (trname_a(k).ne.'!@@NOT@@')
$ call hinterpo
(trr(1,1,k),nid,njd,trn((k-1)*nga*nka+1),nis,njs,
$ nka,idx,idy,cxa,cxb,cxc,cxd,cya,cyb,cyc,cyd,Lam_hint_S)
end do
* Horizontal interpolation (xpus,ypqs) ===> (xpud,ypqd)
call grid_to_grid_coef (xpud,nid,xpus,nis,idu,cua,cub,cuc,cud,
$ Lam_hint_S)
call hinterpo
(uur,nid,njd,uun,nis,njs,nka,
$ idu,idy,cua,cub,cuc,cud,cya,cyb,cyc,cyd,Lam_hint_S)
if (Vertint_L) then
* Compute p0,tt,gz on U grid from analysis
* Horizontal interpolation (xpqs,ypqs) ===> (xpud,ypqd)
call grid_to_grid_coef (xpud,nid,xpqs,nis,idu,cua,cub,cuc,cud,
$ Lam_hint_S)
call hinterpo
(ana_p0u,nid,njd,wlnph,nis,njs,1,
$ idu,idy,cua,cub,cuc,cud,cya,cyb,cyc,cyd,Lam_hint_S)
call hinterpo
(ttru,nid,njd,ttn,nis,njs,nka,
$ idu,idy,cua,cub,cuc,cud,cya,cyb,cyc,cyd,Lam_hint_S)
call hinterpo
(firu(1,nka),nid,njd, fin( 1,nka),nis,njs,1,
$ idu,idy,cua,cub,cuc,cud,cya,cyb,cyc,cyd,Lam_hint_S)
endif
* Horizontal interpolation (xpqs,ypvs) ===> (xpqd,ypvd)
call grid_to_grid_coef (ypvd,njd,ypvs,njs,idu,cua,cub,cuc,cud,
$ Lam_hint_S)
call hinterpo
(vvr,nid,njd,vvn,nis,njs,nka,
$ idx,idu,cxa,cxb,cxc,cxd,cua,cub,cuc,cud,Lam_hint_S)
if (Vertint_L) then
* Compute p0,tt,gz on V grid from analysis
* Horizontal interpolation (xpqs,ypqs) ===> (xpqd,ypvd)
call grid_to_grid_coef (ypvd,njd,ypqs,njs,idu,cua,cub,cuc,cud,
$ Lam_hint_S)
call hinterpo
(ana_p0v,nid,njd,wlnph,nis,njs,1,
$ idx,idu,cxa,cxb,cxc,cxd,cua,cub,cuc,cud,Lam_hint_S)
call hinterpo
(ttrv,nid,njd,ttn,nis,njs,nka,
$ idx,idu,cxa,cxb,cxc,cxd,cua,cub,cuc,cud,Lam_hint_S)
call hinterpo
(firv(1,nka),nid,njd, fin( 1,nka),nis,njs,1,
$ idx,idu,cxa,cxb,cxc,cxd,cua,cub,cuc,cud,Lam_hint_S)
endif
deallocate (idx,idy,idu,cxa,cxb,cxc,cxd,cua,cub,cuc,cud,
$ cya,cyb,cyc,cyd)
allocate (ps(nid,njd),psu(nid,njd),psv(nid,njd),
$ gz_temp(ngd,max(G_nk,nka)))
allocate ( uu3(n1:d1+n2,n3:d2+n4,lnk),
$ vv3(n1:d1+n2,n3:d2+n4,lnk),
$ psd3(n1:d1+n2,n3:d2+n4,lnk),
$ tt3(n1:d1+n2,n3:d2+n4,lnk),
$ tp3(n1:d1+n2,n3:d2+n4,lnk),
$ td3(n1:d1+n2,n3:d2+n4,lnk),
$ fi3(n1:d1+n2,n3:d2+n4,lnk),
$ qq3(n1:d1+n2,n3:d2+n4,lnk),
$ ss3(n1:d1+n2,n3:d2+n4,1),
$ fip3(n1:d1+n2,n3:d2+n4,lnk),
$ pip3(n1:d1+n2,n3:d2+n4,lnk),
$ tr3 (n1:d1+n2,n3:d2+n4,lnk),
$ ww3(n1:d1+n2,n3:d2+n4,lnk),
$ mu3(n1:d1+n2,n3:d2+n4,lnk) )
allocate (topo_temp(n1:d1+n2,n3:d2+n4),
$ topu_temp(n1:d1+n2,n3:d2+n4),
$ topv_temp(n1:d1+n2,n3:d2+n4) )
if (b1) then
do j=n3,n4
do i=n1,n2
topo_temp(i,j) = topo1(i,j)
topu_temp(i,j) = topu1(i,j)
topv_temp(i,j) = topv1(i,j)
enddo
enddo
endif
if (b2) then
do j=n3,n4
do i=n1,n2
topo_temp(i+d1,j+d2)=topo2(i,j)
topu_temp(i+d1,j+d2)=topu2(i,j)
topv_temp(i+d1,j+d2)=topv2(i,j)
enddo
enddo
endif
if (Vertint_L) then
* Preparation for vertical interpolation, grid phi
* Compute vt from analysis using TT and HU
* Compute hydrostatic GZ on NKA analysis levels
do i=1,ngd
gz_temp(i,nka) = fir(i,nka)
enddo
call p0vt2gz_hyb
( gz_temp, ana_pia, ana_pibb, ana_p0,
$ ttr,ngd, nka,.false.,.false.)
* Finally, compute pressure on model surface (PS) using p0 and GZ of analysis
call getp0
( ps, topo_temp,ana_pia,ana_pibb, ana_p0, gz_temp,
$ ttr, ngd, nka,.false.)
* Interpolate VT
call vte_hyb2hyb
(tt3, Geomg_pia, Geomg_pibb,ps,G_nk, ttr,
$ ana_pia,ana_pibb, ana_p0, nka, ngd, 'VT',.false.)
*
* Compute hydrostatic GZ on G_nk model levels
do j=1,njd
do i=1,nid
fi3(i+n1-1,j+n3-1,lnk)=topo_temp(i+n1-1,j+n3-1)
enddo
enddo
call p0vt2gz_hyb
( fi3, Geomg_pia, Geomg_pibb, ps,
$ tt3, ngd, G_nk,.false.,.false. )
*
* Preparation for vertical interpolation, grid U
* Compute hydrostatic GZ on analysis levels
do i=1,ngd
gz_temp(i,nka) = firu(i,nka)
enddo
call p0vt2gz_hyb
( gz_temp, ana_pia, ana_pibb, ana_p0u,
$ ttru,ngd, nka,.false.,.false.)
* Finally, compute pressure on model surface (PS)
call getp0
( psu, topu_temp,ana_pia,ana_pibb, ana_p0,
$ gz_temp, ttru, ngd, nka,.false.)
*
* Preparation for vertical interpolation, grid V
* Compute hydrostatic GZ on analysis levels
do i=1,ngd
gz_temp(i,nka) = firv(i,nka)
enddo
call p0vt2gz_hyb
( gz_temp, ana_pia, ana_pibb, ana_p0v,
$ ttru,ngd, nka,.false.,.false.)
* Finally, compute pressure on model surface (PS)
* psv=ana_p0v
call getp0
( psv, topv_temp,ana_pia,ana_pibb, ana_p0,
$ gz_temp, ttrv, ngd, nka,.false.)
*
* Interpolate UT1
call vte_hyb2hyb
(uu3, Geomg_pia, Geomg_pibb,psu,G_nk, uur,
$ ana_pia,ana_pibb, ana_p0u, nka, ngd, 'UU',.false.)
* Interpolate VT1
call vte_hyb2hyb
(vv3, Geomg_pia, Geomg_pibb,psv,G_nk, vvr,
$ ana_pia,ana_pibb, ana_p0v, nka, ngd, 'VV',.false.)
do j=1,njd
do i=1,nid
qq3(i+n1-1,j+n3-1,lnk)=alog(ps(i,j))
qq3(i+n1-1,j+n3-1, 1)=alog(Pres_ptop*100.)
enddo
enddo
* BCS_predat does not calculate wind related fields. These calculations
* would be done under the space which contains halos
call bcs_predat
(ss3,qq3,pip3,fi3,fip3,tt3,tp3,nid,njd,lnk)
else
* NO VERTINT
do k=1,lnk
do j=1,njd
do i=1,nid
ii=(j-1)*nid+i
fi3(i+n1-1,j+n3-1,k)=fir(ii,k)
tt3(i+n1-1,j+n3-1,k)=ttr(ii,k)
uu3(i+n1-1,j+n3-1,k)=uur(ii,k)
vv3(i+n1-1,j+n3-1,k)=vvr(ii,k)
psd3(i+n1-1,j+n3-1,k)=psdr(ii,k)
td3(i+n1-1,j+n3-1,k)=tdr(ii,k)
tp3(i+n1-1,j+n3-1,k)=tpr(ii,k)
fip3(i+n1-1,j+n3-1,k)=fipr(ii,k)
pip3(i+n1-1,j+n3-1,k)=pipr(ii,k)
qq3(i+n1-1,j+n3-1,k)=qqr(ii,k)
enddo
enddo
enddo
do j=1,njd
do i=1,nid
ii=(j-1)*nid+i
ss3(i+n1-1,j+n3-1,1)=ssr(ii,1)
enddo
enddo
If (.not.Schm_hydro_L.and.nvar.gt.11) then
do k=1,lnk
do j=1,njd
do i=1,nid
ii=(j-1)*nid+i
ww3(i+n1-1,j+n3-1,k)=wwr(ii,k)
mu3(i+n1-1,j+n3-1,k)=mur(ii,k)
enddo
enddo
enddo
endif
endif
*
* Interpolate Tracers and place in bcs space
do 100 n=1,Tr3d_ntr
id = (n-1)*bcs_sz
jj=-1
do k=1,ntra
if (Tr3d_name_S(n).eq.trname_a(k)(1:4)) jj=k
end do
* If data found for this tracer
if (jj.gt.0) then
if (Vertint_L) then
call vte_hyb2hyb
(tr3,Geomg_pia,Geomg_pibb,ps,G_nk,trr(1,1,jj),
$ ana_pia,ana_pibb, ana_p0, nka,ngd,trname_a(jj)(1:2),.false.)
else
do k=1,lnk
do j=1,njd
do i=1,nid
ii=(j-1)*nid+i
tr3(i+n1-1,j+n3-1,k)=trr(ii,k,jj)
enddo
enddo
enddo
endif
* ALWAYS clip tracers to zero after vertical interpolation (Desgagne)
do k=1,lnk
do j=1,njd
do i=1,nid
tr3(i+n1-1,j+n3-1,k)=max(tr3(i+n1-1,j+n3-1,k),0.0)
enddo
enddo
enddo
else
* No data found for this tracer, set to user-defined value
do k=1,lnk
do j=1,njd
do i=1,nid
tr3(i+n1-1,j+n3-1,k)=tr3d_sval(n)
enddo
enddo
enddo
endif
if (.not.Schm_moist_L) then
jjj=-1
* See if it is a humid tracer
do kk = 1,h2o_ntr
if (trname_a(jj)(1:2).eq.h2o_name_S(kk)(1:2)) jjj=kk
enddo
if (jjj.gt.0) then
* If no moist scheme, put humid tracers to zero
do k=1,lnk
do j=1,njd
do i=1,nid
tr3(i+n1-1,j+n3-1,k)=0.0
enddo
enddo
enddo
endif
endif
* Fill the BCS area
if (b1) then
cnt=0
do k=1,lnk
do j=n3,n4
do i=n1,n2
cnt=cnt+1
tr1(cnt + id+ib1) = tr3(i,j,k)
end do
end do
end do
endif
if (b2) then
cnt=0
do k=1,lnk
do j=n3,n4
do i=n1,n2
cnt=cnt+1
tr1(cnt + id+ib2) = tr3(i+d1,j+d2,k)
end do
end do
end do
endif
100 continue
* Place other fields in BCS space
if (b1) then
do k=1,lnk
do j=n3,n4
do i=n1,n2
uu1(i,j,k) = uu3(i,j,k)
vv1(i,j,k) = vv3(i,j,k)
tt1(i,j,k) = tt3(i,j,k)
tp1(i,j,k) = tp3(i,j,k)
fi1(i,j,k) = fi3(i,j,k)
qq1(i,j,k) = qq3(i,j,k)
fip1(i,j,k) = fip3(i,j,k)
pip1(i,j,k) = pip3(i,j,k)
enddo
enddo
enddo
do j=n3,n4
do i=n1,n2
ss1(i,j,1) = ss3(i,j,1)
enddo
enddo
endif
*
if (b2) then
do k=1,lnk
do j=n3,n4
do i=n1,n2
uu2(i,j,k) = uu3(i+d1,j+d2,k)
vv2(i,j,k) = vv3(i+d1,j+d2,k)
tt2(i,j,k) = tt3(i+d1,j+d2,k)
tp2(i,j,k) = tp3(i+d1,j+d2,k)
fi2(i,j,k) = fi3(i+d1,j+d2,k)
qq2(i,j,k) = qq3(i+d1,j+d2,k)
fip2(i,j,k) = fip3(i+d1,j+d2,k)
pip2(i,j,k) = pip3(i+d1,j+d2,k)
end do
end do
end do
do j=n3,n4
do i=n1,n2
ss2(i,j,1) = ss3(i+d1,j+d2,1)
enddo
enddo
endif
if (.not. Vertint_L) then
if (b1) then
do k=1,lnk
do j=n3,n4
do i=n1,n2
psd1(i,j,k) = psd3(i,j,k)
td1(i,j,k) = td3(i,j,k)
enddo
enddo
enddo
endif
if (b2) then
do k=1,lnk
do j=n3,n4
do i=n1,n2
psd2(i,j,k) = psd3(i+d1,j+d2,k)
td2(i,j,k) = td3(i+d1,j+d2,k)
enddo
enddo
enddo
endif
endif
*
*
If (.not.Schm_hydro_L.and.nvar.gt.11) then
if (b1) then
do k=1,lnk
do j=n3,n4
do i=n1,n2
ww1(i,j,k) = ww3(i,j,k)
mu1(i,j,k) = mu3(i,j,k)
enddo
enddo
enddo
endif
if (b2) then
do k=1,lnk
do j=n3,n4
do i=n1,n2
ww2(i,j,k) = ww3(i+d1,j+d2,k)
mu2(i,j,k) = mu3(i+d1,j+d2,k)
enddo
enddo
enddo
endif
endif
deallocate (uur,vvr,psdr,ttr,tpr,tdr,fir,
$ qqr,ssr,fipr,pipr,trr,wwr,mur,
$ wlnph,ps,psu,psv,
$ ana_p0,ana_p0u,ana_p0v,gz_temp,
$ firu,firv,ttru,ttrv)
deallocate (uu3,vv3,psd3,tt3,tp3,td3,fi3,
$ qq3,ss3,fip3,pip3,tr3,ww3,mu3)
deallocate (topo_temp,topu_temp,topv_temp)
*
*-----------------------------------------------------------------------
return
end
*