!-------------------------------------- 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 v4d_cstjok - Interpolate model states at obs. locations and
* output the result to be read by 3D-Var
*
#include "model_macros_f.h"
*
subroutine v4d_cstjok 1,27
*
use v4d_prof
use v4dz
use v4d_interint0
*
* ----------------------------------------------------------------
interface
subroutine v4d_putfld (cdvar,kstatus,gomvar_8,kdim,kobs,indexvar)
*
#include "model_macros_f.h"
*
use v4d_prof
*
implicit none
*
character*2,INTENT(IN) :: cdvar
integer, INTENT(INOUT):: kstatus
integer, INTENT(IN) :: kdim,kobs,indexvar
real*8, pointer, dimension(:,:) :: gomvar_8
*
end subroutine v4d_putfld
end interface
* ----------------------------------------------------------------
*
*author N. Ek
*
*revision
* v3_00 - N.Ek - initial MPI version
* v3_00 M.Tanguay - adapt to Simon's exchange
* v3_01 - M.Tanguay - correction for empty processors
* v3_02 - M.Tanguay - locate HU in tracers
* v3_11 - M.Tanguay - correct relationship between stepob and timestep
* - Add option for profiles done on U-V grids for winds
* v3_30 - Fillion/Tanguay - Adapt diagnostics for LAM
* v3_31 - Tanguay M. - Add OPENMP directives
*
*object
* -----------------------------------------------------------------------
* Each processor do the following at the current bin:
* 1) Conversion from GEM units to 3D-Var units and Reverse Staggering
* 2) Interpolate model profiles to obs. locations
* 3) Write model profiles to be read by 3D-Var
* -----------------------------------------------------------------------
*
*arguments
*
*implicits
#include "glb_ld.cdk"
#include "lctl.cdk"
#include "lun.cdk"
#include "v4dg.cdk"
#include "vt1.cdk"
#include "vt1m.cdk"
#include "ptopo.cdk"
#include "tr3d.cdk"
#include "path.cdk"
#include <clib_interface.cdk>
#include <prof_f.h>
*
* Local variables
* ---------------
integer vmmlod, vmmget, vmmuld, prof_wrrec
external vmmlod, vmmget, vmmuld, prof_wrrec
*
integer nvar, kstatus, pnlkey1(9), err, pnerr, istat, i, j,
% k, n, i1, i2, j1, j2, nbin, npr, inn
*
integer key1(Tr3d_ntr), key1_, key1m(Tr3d_ntr), key1m_
real hut1, hut1m
pointer (pahu1, hut1(LDIST_SHAPE,*)), (pahu1m, hut1m(LDIST_SHAPE,*))
*
real, pointer, dimension(:,:) :: profx,profy,prof2d,prbid
real*8, pointer, dimension(:,:) :: gomu_8,gomv_8,gomt_8,gomq_8,gomps_8
*
real wijk1(LDIST_SHAPE,l_nk),wijk2(LDIST_SHAPE,l_nk),wbid(LDIST_SHAPE,l_nk),
% wijk3(LDIST_SHAPE,l_nk),wijk4(LDIST_SHAPE,l_nk),wij5(LDIST_SHAPE)
*
logical plpr_L
*
character(len=2) :: cljx_S, cljy_S
*
real*8, parameter :: ZERO_8 = 0.0
* ______________________________________________________
*
write(Lun_out,1000) Lctl_step
*
* Nullify pointers
* ----------------
nullify (gomu_8,gomv_8,gomt_8,gomq_8,gomps_8)
*
* Flags for diagnostics
* ---------------------
plpr_L = .false.
plpr_L = plpr_L.and.Lun_out.gt.0
*
* Open dwyf PROF file to write the model profiles
* -----------------------------------------------
if(.not.Pr_wopen_L) then
*
write(cljx_S,'(i2.2)') Ptopo_mycol
write(cljy_S,'(i2.2)') Ptopo_myrow
*
Pr_type2file_S = trim(Path_xchg_S)//'/dwyf_'
% //cljx_S//'_'//cljy_S//'.prof'
*
write(Lun_out,*) 'Opening dwyf MODEL-PROFILE output file'
*
Pr_ihdlout = prof_open (Pr_type2file_S,'WRITE','FILE')
*
if(Pr_ihdlout.le.0) then
write(Lun_out,*) 'Cannot open MODEL-PROFILE output file !'
kstatus = - 99
endif
*
Pr_wopen_L=.true.
*
endif
*
* Recall the dimensions of the fields presented to the interpolation
* ------------------------------------------------------------------
i1=V4dz_i1
i2=V4dz_i2
j1=V4dz_j1
j2=V4dz_j2
*
* Establish at which bin we are
* -----------------------------
nbin = (Lctl_step - Pr_ibin0) / V4dg_stepob + 1
*
* Get fields in memory
* --------------------
pnlkey1(1)= VMM_KEY(ut1)
pnlkey1(2)= VMM_KEY(vt1)
pnlkey1(3)= VMM_KEY(tpt1)
pnlkey1(4)= VMM_KEY(st1)
nvar= 4
*
if(V4dg_tl_L) then
pnlkey1(nvar+1)= VMM_KEY(tpt1m)
pnlkey1(nvar+2)= VMM_KEY(st1m)
nvar= nvar+2
endif
*
err = vmmlod(pnlkey1,nvar)
*
err = VMM_GET_VAR(ut1)
err = VMM_GET_VAR(vt1)
err = VMM_GET_VAR(tpt1)
err = VMM_GET_VAR(st1)
if(V4dg_tl_L) then
err = VMM_GET_VAR(tpt1m)
err = VMM_GET_VAR(st1m)
endif
*
* Load humidity field
* -------------------
key1_ = VMM_KEY (trt1)
do n=1,Tr3d_ntr
key1(n) = key1_ + n
end do
err = vmmlod(key1,Tr3d_ntr)
do n=1,Tr3d_ntr
if (Tr3d_name_S(n).eq.'HU') err = vmmget(key1(n),pahu1,hut1)
end do
*
if(V4dg_tl_L) then
*
* Load TRAJ humidity field
* ------------------------
key1m_ = VMM_KEY (trt1m)
do n=1,Tr3d_ntr
key1m(n) = key1m_ + n
end do
err = vmmlod(key1m,Tr3d_ntr)
do n=1,Tr3d_ntr
if (Tr3d_name_S(n).eq.'HU') err = vmmget(key1m(n),pahu1m,hut1m)
end do
*
end if
*
* Transfer fields
* ---------------
!$omp parallel do
do k=1,l_nk
do j=1,l_nj
do i=1,l_niu
wijk1(i,j,k) = ut1 (i,j,k)
end do
do i=1,l_ni
wijk3(i,j,k) = tpt1(i,j,k)
wijk4(i,j,k) = hut1(i,j,k)
end do
end do
do j=1,l_njv
do i=1,l_ni
wijk2(i,j,k) = vt1 (i,j,k)
end do
end do
end do
!$omp end parallel do
*
!$omp parallel do
do j=1,l_nj
do i=1,l_ni
wij5(i,j) = st1(i,j)
end do
end do
!$omp end parallel do
*
if(plpr_L) then
inn= 0
if (G_lam) then
inn=1
endif
if(Ptopo_myproc.eq.0) write(Lun_out,*) 'BEFORE VARCONV'
call glbstat
(wijk1,'UU',LDIST_DIM,G_nk,1,G_ni-inn,1,G_nj,1,G_nk)
call glbstat
(wijk2,'VV',LDIST_DIM,G_nk,1,G_ni,1,G_nj-1,1,G_nk)
call glbstat
(wijk3,'TP',LDIST_DIM,G_nk,1,G_ni,1,G_nj,1,G_nk)
call glbstat
(wij5 ,'4S',LDIST_DIM, 1,1,G_ni,1,G_nj,1, 1)
call glbstat
(wijk4,'HU',LDIST_DIM,G_nk,1,G_ni,1,G_nj,1,G_nk)
if(Ptopo_myproc.eq.0) write(Lun_out,*) '-----------------------'
endif
*
* --------------------------------------------------------------------------------------
* Conversion from GEM units to 3D-Var units and Reverse Staggering (if .not.V4dg_pruv_L)
* --------------------------------------------------------------------------------------
* Direct (nonlinear)
* ------------------
if(V4dg_di_L) then
*
call v4d_varconv
(wijk1,wijk2,wijk3,wijk4,wij5,LDIST_DIM,l_nk,.FALSE.)
*
* TLM
* ---
elseif(V4dg_tl_L) then
*
call v4d_varconv_tl
(wijk1,wijk2,wijk3,wijk4,wij5,
% tpt1m,hut1m,st1m,LDIST_DIM,l_nk,.FALSE.)
*
end if
*
if(plpr_L) then
if(Ptopo_myproc.eq.0) write(Lun_out,*) 'BEFORE PROFILE'
call glbstat
(wijk1,'UU',LDIST_DIM,G_nk,1,G_ni,1,G_nj,1,G_nk)
call glbstat
(wijk2,'VV',LDIST_DIM,G_nk,1,G_ni,1,G_nj-1,1,G_nk)
call glbstat
(wijk3,'TP',LDIST_DIM,G_nk,1,G_ni,1,G_nj,1,G_nk)
call glbstat
(wij5 ,'4S',LDIST_DIM, 1,1,G_ni,1,G_nj,1, 1)
call glbstat
(wijk4,'HU',LDIST_DIM,G_nk,1,G_ni,1,G_nj,1,G_nk)
if(Ptopo_myproc.eq.0) write(Lun_out,*) '-----------------------'
endif
*
* -------------------------------------------
* Evaluate profiles at observations locations
* -------------------------------------------
*
* --------------------------------
* Contribution from U-V components
* --------------------------------
npr = Pr_l_mv(V3D_UTRU,nbin) % nprof
write(Lun_out,*) 'Evaluate profiles UV at BIN = ',nbin,
% 'Number of profiles = ',npr
*
if(npr.ne.0) allocate(gomu_8(l_nk,npr),gomv_8(l_nk,npr),profx(l_nk,npr),profy(l_nk,npr))
if(npr.ne.0.and.V4dg_pruv_L) allocate(prbid(l_nk,npr))
*
* Interpolation to observation locations using EZSCINT
* ----------------------------------------------------
if(.not.V4dg_pruv_L) then
*
call v4d_uvint0
(profx,profy,Pr_l_mv(V3D_UTRU,nbin)%px,Pr_l_mv(V3D_UTRU,nbin)%py,npr,
% wijk1,wijk2,V4dz_ax,V4dz_ay,V4dz_cx,V4dz_cy,V4dz_wx_8,V4dz_cox_8,V4dz_six_8,V4dz_siy_8,
% i1,i2,j1,j2,l_nk,V4dz_grtypi,V4dz_degree,'UV')
*
else
*
!$omp parallel do
do k=1,l_nk
do j=l_miny,l_maxy
do i=l_minx,l_maxx
wbid(i,j,k) = ZERO_8
end do
end do
end do
!$omp end parallel do
*
call v4d_uvint0
(profx,prbid,Pr_l_mv(V3D_UTRU,nbin)%pxu,Pr_l_mv(V3D_UTRU,nbin)%py,npr,
% wijk1,wbid,V4dz_axu,V4dz_ay,V4dz_cxu,V4dz_cy,V4dz_wxu_8,V4dz_coxu_8,V4dz_sixu_8,V4dz_siy_8,
% i1,i2,j1,j2,l_nk,'U',V4dz_degree,'UV')
*
!$omp parallel do
do k=1,l_nk
do j=l_miny,l_maxy
do i=l_minx,l_maxx
wbid(i,j,k) = ZERO_8
end do
end do
end do
!$omp end parallel do
*
call v4d_uvint0
(prbid,profy,Pr_l_mv(V3D_UTRU,nbin)%px,Pr_l_mv(V3D_UTRU,nbin)%pyv,npr,
% wbid,wijk2,V4dz_ax,V4dz_ayv,V4dz_cx,V4dz_cyv,V4dz_wx_8,V4dz_cox_8,V4dz_six_8,V4dz_siyv_8,
% i1,i2,j1,j2,l_nk,'V',V4dz_degree,'UV')
*
endif
*
if(npr.ne.0) then
*
* Store UU profiles
* -----------------
!$omp parallel do
do j = 1, npr
do k = 1, l_nk
gomu_8(k,j) = dble(profx(k,j))
enddo
enddo
!$omp end parallel do
*
endif
*
call v4d_putfld
('UU',kstatus,gomu_8,l_nk,npr,V3D_UTRU)
*
if(npr.ne.0) deallocate(profx,gomu_8)
*
if(npr.ne.0) then
*
* Store VV profiles
* -----------------
!$omp parallel do
do j = 1, npr
do k = 1, l_nk
gomv_8(k,j) = dble(profy(k,j))
enddo
enddo
!$omp end parallel do
*
endif
*
call v4d_putfld
('VV',kstatus,gomv_8,l_nk,npr,V3D_VTRU)
*
if(npr.ne.0) deallocate(profy,gomv_8)
*
if(npr.ne.0.and.V4dg_pruv_L) deallocate(prbid)
*
* -----------------------------
* Contribution from temperature
* -----------------------------
npr = Pr_l_mv(V3D_TEMP,nbin) % nprof
write(Lun_out,*) 'Evaluate profiles TT at BIN = ',nbin,
% 'Number of profiles = ',npr
*
if(npr.ne.0) allocate(gomt_8(l_nk,npr),profx(l_nk,npr))
*
* Interpolation to observation locations using EZSCINT
* ----------------------------------------------------
call v4d_scint0
(profx,Pr_l_mv(V3D_TEMP,nbin)%px,Pr_l_mv(V3D_TEMP,nbin)%py,npr,
% wijk3,V4dz_ax,V4dz_ay,V4dz_cx,V4dz_cy,V4dz_wx_8,i1,i2,j1,j2,l_nk,
% V4dz_grtypi,V4dz_degree,'TT')
*
if(npr.ne.0) then
*
* Store Temp. profiles
* --------------------
!$omp parallel do
do j = 1, npr
do k = 1, l_nk
gomt_8(k,j) = dble(profx(k,j))
enddo
enddo
!$omp end parallel do
*
endif
*
call v4d_putfld
('TT',kstatus,gomt_8,l_nk,npr,V3D_TEMP)
*
if(npr.ne.0) deallocate(profx,gomt_8)
*
* --------------------------
* Contribution from humidity
* --------------------------
npr = Pr_l_mv(V3D_SPHU,nbin) % nprof
write(Lun_out,*) 'Evaluate profiles HU at BIN = ',nbin,
% 'Number of profiles = ',npr
*
if(npr.ne.0) allocate(gomq_8(l_nk,npr),profx(l_nk,npr))
*
* Interpolation to observation locations using EZSCINT
* ----------------------------------------------------
call v4d_scint0
(profx,Pr_l_mv(V3D_SPHU,nbin)%px,Pr_l_mv(V3D_SPHU,nbin)%py,npr,
% wijk4,V4dz_ax,V4dz_ay,V4dz_cx,V4dz_cy,V4dz_wx_8,i1,i2,j1,j2,l_nk,
% V4dz_grtypi,V4dz_degree,'HU')
*
if(npr.ne.0) then
*
* Store Humidity profiles
* -----------------------
!$omp parallel do
do j = 1, npr
do k = 1, l_nk
gomq_8(k,j) = dble(profx(k,j))
enddo
enddo
!$omp end parallel do
*
endif
*
call v4d_putfld
('HU',kstatus,gomq_8,l_nk,npr,V3D_SPHU)
*
if(npr.ne.0) deallocate(profx,gomq_8)
*
* ----------------------------------
* Contribution from surface pressure
* ----------------------------------
npr = Pr_l_mv(V2D_PSUR,nbin) % nprof
write(Lun_out,*) 'Evaluate profiles PS at BIN = ',nbin,
% 'Number of profiles = ',npr
*
if(npr.ne.0) allocate (gomps_8(1,npr),prof2d(1,npr))
*
* Interpolation to observation locations using EZSCINT
* ----------------------------------------------------
call v4d_scint0
(prof2d,Pr_l_mv(V2D_PSUR,nbin)%px,Pr_l_mv(V2D_PSUR,nbin)%py,npr,
% wij5,V4dz_ax,V4dz_ay,V4dz_cx,V4dz_cy,V4dz_wx_8,i1,i2,j1,j2,1,
% V4dz_grtypi,V4dz_degree,'4S')
*
if(npr.ne.0) then
*
* Store Psurf. profiles
* ---------------------
!$omp parallel do
do j=1,npr
gomps_8(1,j) = dble(prof2d(1,j))
enddo
!$omp end parallel do
*
endif
*
call v4d_putfld
('SP',kstatus,gomps_8,1,npr,V2D_PSUR)
*
if(npr.ne.0) deallocate(prof2d,gomps_8)
*
pnerr = vmmuld(-1,0)
*
write(Lun_out,1001) Lctl_step
call flush(Lun_out)
*
return
*
1000 format(/,'V4D_CSTJOK_TL: Beginning for TIMESTEP = ',I8,
+ /,'=====================================')
1001 format(/,'V4D_CSTJOK_TL: Ending for TIMESTEP = ',I8,
+ /,'=====================================')
*
end