!-------------------------------------- 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 casc_3df_dynp_offline - For reading 3DF file where
* 3DF02 files were written by GEMENT
* Routine to give initial conditions for the model run
*
#include "model_macros_f.h"
*
subroutine casc_3df_dynp_offline (dimgx,dimgy,unf) 1,14
implicit none
*
integer dimgx,dimgy,unf
*
*author
* V. Lee July 2008 (GEM casc_3df_dynp)
*
*revision
* v3_31 - Lee V. - initial version for GEMDM
*
#include "glb_ld.cdk"
#include "dcst.cdk"
#include "geomg.cdk"
#include "ifd.cdk"
#include "ind.cdk"
#include "pres.cdk"
#include "lam.cdk"
#include "ptopo.cdk"
#include "schm.cdk"
#include "tr3d.cdk"
#include "vt1.cdk"
#include "lun.cdk"
#include "itf_phy_buses.cdk"
#include "lctl.cdk"
#include "itf_phy_busind.cdk"
#include "vtopo.cdk"
#include "path.cdk"
*
integer vmmlod,vmmget,vmmuld,longueur,sid3df
external vmmlod,vmmget,vmmuld,longueur,sid3df
*
character*2 md
character*4 nomvar
character*8 dynophy
character*8, dimension (:), pointer :: trname_a
character*15 datev
character*256 fn
logical done,dyn_done
logical dyn_init
integer*8 pnt_trp(Tr3d_ntr)
integer i,j,k,jj,jjj,kk,nia,nja,nka,ntra,ntra1,ni1,nj1,n,err,
$ errop,ofi,ofj,l_in,l_jn,mode,nvar,offg,
$ cnt,errdyn,cumerr,gid,wowp
integer i0,in,j0,jn,ng,keyp_,keyp(Tr3d_ntr)
real topo_temp(l_ni,l_nj)
integer, dimension (: ), pointer :: idx,idy,nks
real psmin, psmax, psmin_glob, psmax_glob
real*8, dimension (: ), pointer ::
$ xpaq,ypaq,xpau,ypav,
$ cxa,cxb,cxc,cxd,cya,cyb,cyc,cyd
real, dimension (:,:), pointer ::
$ uun,vvn,ttn,ps
real, dimension (:,:,:), pointer ::
$ uur,vvr,ttr,trn,trr
real trp
pointer (patrp, trp(LDIST_SHAPE,*))
real*8 xpxext(0:dimgx+1), ypxext(0:dimgy+1)
*-----------------------------------------------------------------------
*
if (Lun_debug_L) write (Lun_out,1000)
*
keyp_ = VMM_KEY (trt1)
do k=1,Tr3d_ntr
keyp(k) = keyp_ + k
end do
err = vmmlod(keyp,Tr3d_ntr)
do k=1,Tr3d_ntr
err = vmmget(keyp(k),patrp,trp)
pnt_trp(k) = patrp
end do
*
* Positional parameters on extended global grid
*
do i=1,dimgx
xpxext(i) = G_xg_8(i)
end do
xpxext(0) = xpxext(1) - (xpxext(2)-xpxext(1))
xpxext(dimgx+1) = xpxext(dimgx) + (xpxext(dimgx)-xpxext(dimgx-1))
*
do i=1,dimgy
ypxext(i) = G_yg_8(i)
end do
ypxext(0) = ypxext(1) - (ypxext(2)-ypxext(1))
ypxext(dimgy+1) = ypxext(dimgy) + (ypxext(dimgy)-ypxext(dimgy-1))
*
* Read all needed files and construct the source domain for
* the horizontal interpolation
*
nia = ifd_niaf - ifd_niad + 1
nja = ifd_njaf - ifd_njad + 1
nullify(xpaq,xpau,ypaq,ypav,trname_a)
nullify(uun,vvn,ttn,ps)
nullify(uur,vvr,ttr,trn,trr)
*
if (associated(xpaq)) deallocate(xpaq)
if (associated(ypaq)) deallocate(ypaq)
if (associated(xpau)) deallocate(xpau)
if (associated(ypav)) deallocate(ypav)
allocate (xpaq(nia), ypaq(nja), xpau(nia), ypav(nja))
*
datev= Lam_runstrt_S
*
ntra = 0
err = 0
*
* wowp = 2 ===> input data has seen the physics
* wowp = 1 ===> input data just after dynamics (no physics)
* We prefer to initialize uup, vvp etc... with wowp=2 status.
*
wowp = 3
48 wowp = wowp - 1
if (wowp.lt.1) then
write (6,204)
err = -1
goto 999
endif
*
write (md,'(i2.2)') wowp
done = .false.
dyn_init = .false.
*
do n=1,ifd_nf
*
ofi = ifd_minx(n)-1
ofj = ifd_miny(n)-1
if (ifd_needit(n)) then
*
errdyn = -1
dyn_done = .false.
*
fn = trim(Path_ind_S)//'/3df'//md//
$ '_'//datev//'_'//ifd_fnext(n)
open (unf,file=fn(1:longueur(fn)),access='SEQUENTIAL',
$ form='UNFORMATTED',status='OLD',iostat=errop)
if (Lun_debug_L)
$ write(Lun_out,*) 'Opening',fn(1:longueur(fn)),'err=',errop
if (errop.ne.0) goto 33
*
* Use first file to establish 3D grid dimensions and geo-references
* of all input staggered grids (xpaq, ypaq, xpau and ypva).
*
55 if (dyn_done) goto 33
err = sid3df
(xpaq,ypaq,xpau,ypav,unf,done,nia,nja,
$ nka,nvar,ntra1)
if (err.lt.0) goto 33
*
read (unf,end=1010,err=1010) dynophy,cnt,mode
*
if (dynophy.eq.'DYNAMICS') then
*
ntra=ntra1
if (.not.dyn_init) then
if (associated(uun)) deallocate(uun)
if (associated(vvn)) deallocate(vvn)
if (associated(ttn)) deallocate(ttn)
if (associated(trn)) deallocate(trn)
if (associated(trname_a)) deallocate(trname_a)
allocate (
$ uun(nia*nja,nka), vvn(nia*nja,nka),
$ ttn(nia*nja,nka),
$ trn(nia*nja,nka,ntra), trname_a(ntra) )
dyn_init = .true.
endif
*
cumerr=0
call filmup
( ttn,ifd_niad,ifd_niaf,ifd_njad,
$ ifd_njaf, nka,unf,ofi,ofj,cumerr )
if (ntra.gt.0) then
call filuptr
( trn,ifd_niad,ifd_niaf,ifd_njad,ifd_njaf,
$ nka,unf,ofi,ofj,Tr3d_name_S,Tr3d_ntr,
$ trname_a,ntra,cumerr )
endif
call filmup
( uun,ifd_niad,ifd_niaf,ifd_njad,ifd_njaf,
$ nka,unf,ofi,ofj,cumerr )
call filmup
( vvn,ifd_niad,ifd_niaf,ifd_njad,ifd_njaf,
$ nka,unf,ofi,ofj,cumerr )
errdyn = cumerr
dyn_done = .true.
else
* Unrecognizable marker
write (6,205) dynophy
goto 1010
endif
*
33 close (unf)
*
if ((errdyn.lt.0).and.(wowp.gt.1)) goto 48
*
err = err + errdyn
done = .true.
if (err.lt.0) then
write (6,203) fn(1:longueur(fn)),Ptopo_myproc
goto 999
endif
endif
end do
*
999 call gem_stop
('casc_3df_dynp_offline',err)
*
* Copy target topography field from geofld (see geodata.ftn) unless
* "growing" topography is used, in which case the current model
* topography is retained.
*
if (Vtopo_L .and. Lctl_step > Vtopo_start) then
topo_temp = Ind_topo(1:l_ni,1:l_nj)
else
do gid=1,P_bgeo_top
if (geonm(gid,1).eq.'MF') then
offg = geopar(gid,1)
cnt = 0
do j=1,l_nj
do i=1,l_ni
cnt=cnt+1
topo_temp(i,j)=dble(geofld(offg +cnt-1))*Dcst_grav_8
enddo
enddo
endif
enddo
endif
*
* Establish geo-references of model target horizontal grids
* (xp1, yp1).
i0 = 1
j0 = 1
in = l_ni
jn = l_nj
ni1 = in - i0 + 1
nj1 = jn - j0 + 1
*
if (associated(uur)) deallocate(uur)
if (associated(vvr)) deallocate(vvr)
if (associated(ttr)) deallocate(ttr)
if (associated(trr)) deallocate(trr)
allocate ( uur(ni1,nj1,nka), vvr(ni1,nj1,nka),
$ ttr(ni1,nj1,nka),
$ trr(ni1*nj1,nka,ntra))
*
ofi = l_i0 - 1
ofj = l_j0 - 1
*
* Horizontal interpolation (xpaq,ypaq) ===> (xp1,yp1) PHI GRID
*
if (associated(idx)) deallocate(idx)
if (associated(idy)) deallocate(idy)
if (associated(cxa)) deallocate(cxa)
if (associated(cxb)) deallocate(cxb)
if (associated(cxc)) deallocate(cxc)
if (associated(cxd)) deallocate(cxd)
if (associated(cya)) deallocate(cya)
if (associated(cyb)) deallocate(cyb)
if (associated(cyc)) deallocate(cyc)
if (associated(cyd)) deallocate(cyd)
allocate (idx(l_ni), idy(l_nj))
allocate (cxa(l_ni),cxb(l_ni),cxc(l_ni),cxd(l_ni),
$ cya(l_nj),cyb(l_nj),cyc(l_nj),cyd(l_nj))
*
call grid_to_grid_coef (xpxext(l_i0),ni1,
$ xpaq,nia,idx,cxa,cxb,cxc,cxd,Lam_hint_S)
call grid_to_grid_coef (ypxext(l_j0),nj1,
$ ypaq,nja,idy,cya,cyb,cyc,cyd,Lam_hint_S)
*
* Interpolate topography from input GZ
* Assume topography is not the same as analysis, always interpolate
*
call hinterpo
( ttr,ni1,nj1, ttn,nia,nja,nka,
$ idx,idy,cxa,cxb,cxc,cxd,cya,cyb,cyc,cyd,Lam_hint_S)
*
* Compute p0, surface pressure from analysis
*
if (associated(ps)) deallocate (ps)
allocate (ps(ni1,nj1))
*
* Humidity is in first cube of trr(1,1,1)
*
do kk=1,ntra
if (trname_a(kk).ne.'!@@NOT@@') then
call hinterpo
(trr(1,1,kk),ni1,nj1,trn(1,1,kk),nia,nja,nka,
$ idx,idy,cxa,cxb,cxc,cxd,cya,cyb,cyc,cyd,Lam_hint_S)
if (trname_a(kk).eq.'P0') then
do j=1,nj1
do i=1,ni1
ps(i,j) = trr(i+(j-1)*ni1,1,kk)
enddo
enddo
endif
endif
end do
*
* Horizontal interpolation (xpaq,ypaq) ===> (xp1,yp1) PHI GRID
*
call hinterpo
(uur,ni1,nj1,uun,nia,nja,nka,
$ idx,idy,cxa,cxb,cxc,cxd,cya,cyb,cyc,cyd,Lam_hint_S)
*
* Horizontal interpolation (xpaq,ypaq) ===> (xp1,yp1) PHI GRID
*
call hinterpo
(vvr,ni1,nj1,vvn,nia,nja,nka,
$ idx,idy,cxa,cxb,cxc,cxd,cya,cyb,cyc,cyd,Lam_hint_S)
*
* Allocate surface pressures for scalar,U,V grid
*
ng = ni1*nj1
*
*
* Compute hydrostatic GZ on model's levels
*
do k=1,G_nk
do j=1,l_nj
do i=1,l_ni
Ind_fi(i,j,k) = 0.0
Ind_t (i,j,k) = ttr(i,j,1)
end do
end do
end do
*
do k=1,G_nk
do j=1,l_nj
do i=1,l_ni
Ind_u(i,j,k) = uur(i,j,1)
Ind_v(i,j,k) = vvr(i,j,1)
end do
end do
end do
*
psmin = ps(1,1)
psmax = ps(1,1)
do j=1,l_nj
do i=1,l_ni
psmin = min( psmin, ps(i,j) )
psmax = max( psmax, ps(i,j) )
Ind_q(i,j,g_nk) = alog(ps(i,j))
Ind_q(i,j,1 ) = alog(Pres_ptop*100.)
enddo
enddo
*
* TRACERS
*
do 200 n=1,Tr3d_ntr
patrp = pnt_trp(n)
jj=-1
* If data found for this tracer
do k=1,ntra
if (Tr3d_name_S(n).eq.trname_a(k)(1:4)) jj=k
end do
if ( jj.gt.0 ) then
* ALWAYS clip tracers to zero after vertical interpolation (Desgagne)
do k=1,G_nk
do j=1,l_nj
do i=1,l_ni
trp(i,j,k) = max(trr(i+(j-1)*l_ni,1,jj),0.0)
end do
end do
end do
else
* No data found for this tracer, set to user-defined value.
do k=1,G_nk
do j=1,l_nj
do i=1,l_ni
trp(i,j,k) = Tr3d_sval(n)
end do
end do
end do
endif
* If no moist scheme, put humid tracers to zero
if (.not.Schm_moist_L) then
jjj=-1
* See if it is a humid tracer
do kk = 1,h2o_ntr
if (Tr3d_name_S(n).eq.h2o_name_S(kk)) jjj=kk
enddo
if (jjj.gt.0) then
do k=1,G_nk
do j=1,l_nj
do i=1,l_ni
trp(i,j,k) = 0.0
enddo
enddo
enddo
endif
endif
200 continue
*
* Copy topography into vmm field
*
do j=1,l_nj
do i=1,l_ni
Ind_topo(i,j) = topo_temp(i,j)
enddo
enddo
*
if (Lun_debug_L) then
write(Lun_out,100)
write(Lun_out,101) datev,wowp
write(Lun_out,100)
endif
*
call rpn_comm_allreduce(psmin,psmin_glob,1,"MPI_REAL","MPI_MIN",
$ "grid",err)
call rpn_comm_allreduce(psmax,psmax_glob,1,"MPI_REAL","MPI_MAX",
$ "grid",err)
psmin=psmin_glob
psmax=psmax_glob
*
if ( Ptopo_myproc.eq.0 ) then
write(6,*)'PSMIN = ',PSMIN,' PSMAX = ',PSMAX,
$ ' PSMINMAX = ',0.5*(PSMIN+PSMAX),' (PASCAL)'
endif
*
Pres_surf = dble(0.5*(psmin+psmax))
Pres_top = dble(Pres_ptop*100.)
*
call v4d_indata3
()
call set_dync
*
call predat
()
if ( .not. Schm_hydro_L ) then
Ind_mul = 0.
Ind_qp = 0.
endif
*
100 format (' ',65('*'))
101 format (' (casc_3df_dynp_offline) JUST READ INIT DATA FOR DATE: ',a15,1x,i3)
203 format (/' PROBLEM WITH FILE: ',a,', PROC#:',i4,' --ABORT--'/)
204 format (/' NO DATA IN casc_3df_dynp_offline --ABORT--'/)
205 format (/' Unrecognizable tag found: ',a,'?'/)
1000 format(
+3X,'READING DATA IN (S/R casc_3df_dynp_offline)')
*
*-----------------------------------------------------------------------
return
1010 write (6,203) fn(1:longueur(fn)),Ptopo_myproc
call gem_stop
('casc_3df_dynp_offline',-1)
return
end