!-------------------------------------- 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 --------------------------------------
!
subroutine getprof 1,7
use mod4dv
, only : maxxy
USE procs_topo
USE obstag
,only : ObsTagLoc
implicit none
*
*author P.Gauthier
*
*revision M.Tanguay
* Simon Pellerin : . Get distributed profiles from TLM
* . Book keeping of PEs, model and observation
* tags and for further communications of
* adjoint profiles (putprofad)
* . Introduction of events (EVN_*) protocol
* S. Pellerin, ARMA, August 2008
* .Added calls to 'tmg_*' subroutines
* Bin He, ARMA, Jan. 2009
* .Implemented MPI to 3/4DVAR.
*
* Correction dlbuff2d
* -------------------
*
#include "comdim.cdk"
#include "comcva.cdk"
#include "comlun.cdk"
#include "comdimo.cdk"
#include "commvo.cdk"
#include "comvfiles.cdk"
include 'prof_f.h'
*
* Local variables
*
integer :: ievent,prof_wrrec,prof_size,ibid_nlev,nprof
character(len=512) :: clprof
character(len=2) :: cljx, cljy
integer jobs, jk, ibprm, ibint, ibreal,ibdbl,jbit,ix,iy, ipes
integer ihdl, prof_rdrec, istat,prof_bitptrn,iobs,ixy,jx,jy
integer, pointer, dimension(:) :: iotag,imtag
real*8, pointer:: dlbuff(:,:), dlbuff2d(:)
nullify(iotag,imtag,dlbuff,dlbuff2d)
*
* Opening Prof_File
*
write(nulout,*) 'Entering GETPROF for simulation ',nsim3d
ix = 0
iy = 0
jx = myidx
jy = myidy
write(cljx,'(i2.2)') jx
write(cljy,'(i2.2)') jy
! Building (jx,jy) bit encoded PEs
ipes = 0
ipes = ishft(jx,16)
call mvbits(jy,0,16,ipes,0)
clprof = trim(CEXC4DV) // '/dwyf_'//cljx//'_'//cljy//'.prof'
call tmg_start(82,'I/O_WAIT')
IHDL = PROF_OPEN(clprof,'READ','FILE')
call tmg_stop(82)
if(ihdl > 0) then
record: do
call tmg_start(81,'PROF_R+W')
ISTAT = prof_rdrec(ihdl)
call tmg_stop(81)
if(istat /= 0) exit record
istat = prof_gvar(ihdl,ievent,prm_evnt)
if(ievent == nsim3d) then
istat = prof_gvar(ihdl,maxxy,prm_pexy) + istat
call mvbits(maxxy,16,16,ix,0)
call mvbits(maxxy,0,16,iy,0)
istat = prof_size(ihdl,ibid_nlev,nprof)
if(nprof.gt.0) then
istat = prof_gvar(ihdl,iotag,v2d_otag) + istat
istat = prof_gvar(ihdl,imtag,v2d_mtag) + istat
do jobs = 1, size(iotag)
iobs = iotag(jobs)
mtag(iobs) = imtag(jobs)
end do
ibprm = 0
ibint = 0
ibreal = 0
ibdbl = 0
istat = prof_bitptrn(ihdl,ibprm,ibint,ibreal,ibdbl)
do jbit = 0,31
if(btest(ibdbl,jbit))
& call getfld
(jbit,iotag,size(iotag))
enddo
endif
else
istat = prof_close(ihdl,.false.)
write(nulout,*)
& 'GETPROF - ERROR - ABNORMAL TERMINATION : ievent = '
& ,ievent
clprof = trim(CEXC4DV) // '/dwya_'//cljx//'_'//cljy/
& /'.prof'
IHDL = PROF_OPEN(clprof,'WRITE','FILE')
ISTAT = PROF_PVAR(IHDL,EVN_FERR,PRM_EVNT)
call tmg_start(81,'PROF_R+W')
istat = prof_wrrec(ihdl)
istat = prof_close(ihdl)
call tmg_stop(81)
call abort3d
(nulout,'GETPROF')
endif
enddo record
istat = prof_close(ihdl,.true.)
else
write(nulout,*)'GETPROF - ERROR: Problem opening file ',clprof
clprof = trim(CEXC4DV) // '/dwya_'//cljx//'_'//cljy/
& /'.prof'
IHDL = PROF_OPEN(clprof,'WRITE','FILE')
ISTAT = PROF_PVAR(IHDL,EVN_FERR,PRM_EVNT)
call tmg_start(81,'PROF_R+W')
istat = prof_wrrec(ihdl)
istat = prof_close(ihdl)
call tmg_stop(81)
call abort3d
(nulout,'GETPROF')
endif
*
write(nulout,*) 'Leaving GETPROF'
call vflush
(nulout)
contains
subroutine getfld(kbit,kotag,kobs) 1
*
integer, intent(in) :: kbit, kobs
integer, intent(in), dimension(kobs) :: kotag
*
* Local variables
integer :: mobs ,iobs
*
select case (kbit)
case(V3D_UTRU)
istat = prof_gvar(IHDL,DLBUFF,V3D_UTRU)
do jobs = 1, kobs
iobs = kotag(jobs)
mobs = ObsTagLoc(iobs)
do jk = 1, nflev
gomu(jk,mobs) = dlbuff(jk,jobs)
end do
end do
case(V3D_VTRU)
istat = prof_gvar(IHDL,DLBUFF,V3D_VTRU)
do jobs = 1, kobs
iobs = kotag(jobs)
mobs = ObsTagLoc(iobs)
do jk = 1, nflev
gomv(jk,mobs) = dlbuff(jk,jobs)
end do
end do
case(V3D_TEMP)
istat = prof_gvar(IHDL,DLBUFF,V3D_TEMP)
do jobs = 1, kobs
iobs = kotag(jobs)
mobs = ObsTagLoc(iobs)
do jk = 1, nflev
gomt(jk,mobs) = dlbuff(jk,jobs)
end do
end do
case(V3D_SPHU)
istat = prof_gvar(IHDL,DLBUFF,V3D_SPHU)
do jobs = 1, kobs
iobs = kotag(jobs)
mobs = ObsTagLoc(iobs)
do jk = 1, nflev
gomq(jk,mobs) = dlbuff(jk,jobs)
end do
end do
case(V2D_PSUR)
istat = prof_gvar(IHDL,DLBUFF2d,V2D_PSUR)
do jobs = 1, kobs
iobs = kotag(jobs)
mobs = ObsTagLoc(iobs)
gomps(1,mobs) = dlbuff2d(jobs)
end do
case(V2D_TGRN)
istat = prof_gvar(IHDL,DLBUFF2d,V2D_TGRN)
do jobs = 1, kobs
iobs = kotag(jobs)
mobs = ObsTagLoc(iobs)
gomtgr(1,mobs) = 0.
end do
end select
*
end subroutine getfld
end subroutine getprof