!-------------------------------------- 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 vspng_drv_tl - TLM of vspng_drv
*
#include "model_macros_f.h"
*
subroutine vspng_drv_tl ( F_u , F_v , F_psd , F_tp , F_ip , 1,36
% F_td , F_t , F_it , F_plt , F_q ,
% F_pip , F_w , F_qp , F_mu , F_mul , F_s ,
% F_um , F_vm, F_psdm, F_tpm , F_ipm ,
% F_tdm , F_tm, F_itm , F_pltm, F_qm ,
% F_pipm, F_wm, F_qpm , F_mum , F_mulm, F_sm,
% F_xfis, DIST_DIM,Nk)
*
implicit none
*
integer DIST_DIM, Nk
real F_u (DIST_SHAPE,Nk), F_v (DIST_SHAPE,Nk),
% F_psd(DIST_SHAPE,Nk), F_tp(DIST_SHAPE,Nk),
% F_ip (DIST_SHAPE,Nk), F_td(DIST_SHAPE,Nk),
% F_t (DIST_SHAPE,Nk), F_it(DIST_SHAPE,Nk),
% F_plt(DIST_SHAPE,Nk), F_q (DIST_SHAPE,Nk),
% F_pip(DIST_SHAPE,Nk), F_w (DIST_SHAPE,Nk),
% F_qp (DIST_SHAPE,Nk), F_mu(DIST_SHAPE,Nk),
% F_mul(DIST_SHAPE,Nk), F_s(DIST_SHAPE), F_xfis(DIST_SHAPE)
*
real F_um (DIST_SHAPE,Nk), F_vm (DIST_SHAPE,Nk),
% F_psdm(DIST_SHAPE,Nk), F_tpm(DIST_SHAPE,Nk),
% F_ipm (DIST_SHAPE,Nk), F_tdm(DIST_SHAPE,Nk),
% F_tm (DIST_SHAPE,Nk), F_itm(DIST_SHAPE,Nk),
% F_pltm(DIST_SHAPE,Nk), F_qm (DIST_SHAPE,Nk),
% F_pipm(DIST_SHAPE,Nk), F_wm (DIST_SHAPE,Nk),
% F_qpm (DIST_SHAPE,Nk), F_mum(DIST_SHAPE,Nk),
% F_mulm(DIST_SHAPE,Nk), F_sm(DIST_SHAPE)
*
*author
* M.Tanguay
*
*revision
* v2_21 - Tanguay M. - initial MPI version
* v2_31 - Tanguay M. - TLM of control for sponge on momentum and
* T', T'lin & hence T on Vspng_nk levels +
* sponge on top level only on all other variables
* v3_00 - Tanguay M. - correction as in vspng_drv
* v3_01 - Laroche/Tanguay - correction nj as in vspng_drv
* v3_01 - Tanguay M. - adapt to Vspng_rwnd_L
* v3_03 - Tanguay M. - Adjoint NoHyd configuration
* v3_11 - Tanguay M. - AIXport+Opti+OpenMP for TLM-ADJ
* v3_20 - Tanguay M. - Adjoint of variable higher order diffusion operator
* - Option of storing instead of redoing TRAJ
* v3_30 - Tanguay M. - adjust TL/AD to Vspng_zmean_L
* v3_31 - Tanguay M. - add Hzd_rwnd_L
*
*object
* see id section
*
*arguments
*
*implicits
#include "glb_ld.cdk"
#include "dcst.cdk"
#include "cstv.cdk"
#include "schm.cdk"
#include "geomg.cdk"
#include "trp.cdk"
#include "vspng.cdk"
#include "hzd.cdk"
#include "opr.cdk"
#include "v4dr.cdk"
*
integer i, j, jj, k, nkspng
real*8 HALF_8,TWO_8,c_8,c1_8
parameter( HALF_8 = 0.5 )
parameter( TWO_8 = 2.0 )
*
real work(LDIST_SHAPE,Nk),work_m(LDIST_SHAPE,Nk)
real tmean(l_nj,Nk)
real*8, dimension (trp_12emax*G_ni*Vspng_nk) :: aix_8,bix_8,cix_8,dix_8
real*8, dimension (trp_22emax*G_nj*Vspng_nk) :: aiy_8,biy_8,ciy_8
real*8 cy_8(l_nj+1), xp0_8(G_ni), yp0_8(G_nj)
* ______________________________________________________
*
if (V4dr_redotr_L) call gem_stop
('STOP in VSPNG_DRV_TL: REDO Not done for Vspng_zmean_L',-1)
*
do i = 1, G_ni
xp0_8 (i) = G_xg_8(i+1) - G_xg_8(i)
end do
do j = 1, G_nj
yp0_8 (j) = sin(G_yg_8(j+1))-sin(G_yg_8(j))
end do
*
* Momentum
* ~~~~~~~~
if ( Cstv_uvdf_8 .gt. 0.0 ) then
*
* Substract the mean for the zonal component if wanted
*
if(Vspng_zmean_L)
$ call vspng_zmean
(F_u,F_u,tmean,DIST_DIM,Nk,.true.)
*
*
* Horizontal Momentum
* ~~~~~~~~~~~~~~~~~~~
do j = 1, l_nj+1
cy_8(j) = G_yg_8(l_j0+j-1)
end do
c_8 = Cstv_uvdf_8*Cstv_dt_8/(Dcst_rayt_8*Dcst_rayt_8)
if (Hzd_difva_L) then
call vspng_abc
( aix_8,bix_8,cix_8,dix_8,aiy_8,biy_8,ciy_8,c_8,
$ cy_8,Hzd_xp0_8,Hzd_xp2_8,Opr_opsyp0_8,Hzd_yp2su_8,
$ G_ni,G_nj,G_nj)
else
call vspng_abc
( aix_8,bix_8,cix_8,dix_8,aiy_8,biy_8,ciy_8,c_8,
$ cy_8,Hzd_xp0_8,Hzd_xp2_8,Opr_opsyp0_8,Opr_opsyp2_8,
$ G_ni,G_nj,G_nj)
endif
*
* TRAJECTORY
* ----------
if (V4dr_redotr_L) then
call vspng_del2
( F_um, xp0_8, Opr_opsyp0_8(G_nj+1),
$ aix_8,bix_8,cix_8,dix_8,aiy_8,biy_8,ciy_8,
$ LDIST_DIM,Vspng_nk,trp_12emax,trp_22emax,G_nj)
endif
*
* TLM
* ---
call vspng_del2
( F_u, xp0_8, Opr_opsyp0_8(G_nj+1),
$ aix_8,bix_8,cix_8,dix_8,aiy_8,biy_8,ciy_8,
$ LDIST_DIM,Vspng_nk,trp_12emax,trp_22emax,G_nj)
*
do j = 1, l_nj+1
jj = l_j0+j-1
cy_8(j) = cos((G_yg_8(jj+1)+G_yg_8(jj)) * HALF_8) **TWO_8
end do
call vspng_abc
( aix_8,bix_8,cix_8,dix_8,aiy_8,biy_8,ciy_8,c_8,
$ cy_8,Opr_opsxp0_8,Opr_opsxp2_8,Hzd_yp0_8,Hzd_yp2_8,
$ G_ni,G_nj,G_njv)
*
* TRAJECTORY
* ----------
if (V4dr_redotr_L) then
call vspng_del2
( F_vm, Opr_opsxp0_8(G_ni+1), yp0_8,
$ aix_8,bix_8,cix_8,dix_8,aiy_8,biy_8,ciy_8,
$ LDIST_DIM,Vspng_nk,trp_12emax,trp_22emax,G_njv)
endif
*
* TLM
* ---
call vspng_del2
( F_v, Opr_opsxp0_8(G_ni+1), yp0_8,
$ aix_8,bix_8,cix_8,dix_8,aiy_8,biy_8,ciy_8,
$ LDIST_DIM,Vspng_nk,trp_12emax,trp_22emax,G_njv)
*
* Add back the mean for the zonal component
*
if (Vspng_zmean_L)
$ call vspng_zmean
(F_u,F_u,tmean,DIST_DIM,Nk,.false.)
*
endif
*
if ( Cstv_phidf_8 .gt. 0.0) then
*
c_8 = Cstv_phidf_8*Cstv_dt_8/(Dcst_rayt_8*Dcst_rayt_8)
do j = 1, l_nj+1
cy_8(j) = G_yg_8(l_j0+j-1)
end do
if (Hzd_difva_L) then
call vspng_abc
(aix_8,bix_8,cix_8,dix_8,aiy_8,biy_8,ciy_8,c_8,
$ cy_8,Opr_opsxp0_8,Opr_opsxp2_8,Opr_opsyp0_8,
$ Hzd_yp2su_8,G_ni,G_nj,G_nj)
else
call vspng_abc
(aix_8,bix_8,cix_8,dix_8,aiy_8,biy_8,ciy_8,c_8,
$ cy_8,Opr_opsxp0_8,Opr_opsxp2_8,Opr_opsyp0_8,
$ Opr_opsyp2_8,G_ni,G_nj,G_nj)
endif
*
* Vertical motion
* ~~~~~~~~~~~~~~~
*
* TRAJECTORY
* ----------
if (V4dr_redotr_L) then
call vspng_del2
(F_psdm,
$ Opr_opsxp0_8(G_ni+1),Opr_opsyp0_8(G_nj+1),
$ aix_8,bix_8,cix_8,dix_8,aiy_8,biy_8,ciy_8,
$ LDIST_DIM,Vspng_nk,trp_12emax,trp_22emax,G_nj)
endif
*
* TLM
* ---
call vspng_del2
(F_psd,
$ Opr_opsxp0_8(G_ni+1),Opr_opsyp0_8(G_nj+1),
$ aix_8,bix_8,cix_8,dix_8,aiy_8,biy_8,ciy_8,
$ LDIST_DIM,Vspng_nk,trp_12emax,trp_22emax,G_nj)
*
* Total divergence
* ~~~~~~~~~~~~~~~~
*
* TRAJECTORY
* ----------
if (V4dr_redotr_L) then
call vspng_del2
(F_tdm,
$ Opr_opsxp0_8(G_ni+1),Opr_opsyp0_8(G_nj+1),
$ aix_8,bix_8,cix_8,dix_8,aiy_8,biy_8,ciy_8,
$ LDIST_DIM,Vspng_nk,trp_12emax,trp_22emax,G_nj)
endif
*
* TLM
* ---
call vspng_del2
(F_td,
$ Opr_opsxp0_8(G_ni+1),Opr_opsyp0_8(G_nj+1),
$ aix_8,bix_8,cix_8,dix_8,aiy_8,biy_8,ciy_8,
$ LDIST_DIM,Vspng_nk,trp_12emax,trp_22emax,G_nj)
*
* The temperature: T', T'lin & hence T
* ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
*
* TRAJECTORY
* ----------
if (V4dr_redotr_L) then
call vspng_del2
(F_tpm,
$ Opr_opsxp0_8(G_ni+1),Opr_opsyp0_8(G_nj+1),
$ aix_8,bix_8,cix_8,dix_8,aiy_8,biy_8,ciy_8,
$ LDIST_DIM,Vspng_nk,trp_12emax,trp_22emax,G_nj)
endif
*
* TLM
* ---
call vspng_del2
(F_tp,
$ Opr_opsxp0_8(G_ni+1),Opr_opsyp0_8(G_nj+1),
$ aix_8,bix_8,cix_8,dix_8,aiy_8,biy_8,ciy_8,
$ LDIST_DIM,Vspng_nk,trp_12emax,trp_22emax,G_nj)
*
* TRAJECTORY
* ----------
if (V4dr_redotr_L) then
call vspng_del2
(F_pltm,
$ Opr_opsxp0_8(G_ni+1),Opr_opsyp0_8(G_nj+1),
$ aix_8,bix_8,cix_8,dix_8,aiy_8,biy_8,ciy_8,
$ LDIST_DIM,Vspng_nk,trp_12emax,trp_22emax,G_nj)
endif
*
* TLM
* ---
call vspng_del2
(F_plt,
$ Opr_opsxp0_8(G_ni+1),Opr_opsyp0_8(G_nj+1),
$ aix_8,bix_8,cix_8,dix_8,aiy_8,biy_8,ciy_8,
$ LDIST_DIM,Vspng_nk,trp_12emax,trp_22emax,G_nj)
*
!$omp parallel
*
if (V4dr_redotr_L) then
!$omp do
do k=1, Vspng_nk
do j=1, l_nj
do i=1, l_ni
*
* TRAJECTORY
* ----------
F_tm(i,j,k) = F_tpm(i,j,k) + Cstv_tstr_8
*
end do
end do
end do
!$omp enddo
endif
*
!$omp do
do k=1, Vspng_nk
do j=1, l_nj
do i=1, l_ni
*
* TLM
* ---
F_t(i,j,k) = F_tp(i,j,k)
*
end do
end do
end do
!$omp enddo
*
!$omp end parallel
*
* Mass related fields
* ~~~~~~~~~~~~~~~~~~~
if (Vspng_uvwdt_L) then
nkspng = Vspng_nk
Vspng_nk = 1
*
if (Hzd_difva_L) then
call vspng_abc
(aix_8,bix_8,cix_8,dix_8,aiy_8,biy_8,
$ ciy_8,c_8,cy_8,Opr_opsxp0_8,Opr_opsxp2_8,
$ Opr_opsyp0_8,Hzd_yp2su_8,G_ni,G_nj,G_nj)
else
call vspng_abc
(aix_8,bix_8,cix_8,dix_8,aiy_8,biy_8,
$ ciy_8,c_8,cy_8,Opr_opsxp0_8,Opr_opsxp2_8,
$ Opr_opsyp0_8,Opr_opsyp2_8,G_ni,G_nj,G_nj)
endif
endif
*
* TRAJECTORY
* ----------
if (V4dr_redotr_L) then
call vspng_del2
(F_ipm,
$ Opr_opsxp0_8(G_ni+1),Opr_opsyp0_8(G_nj+1),
$ aix_8,bix_8,cix_8,dix_8,aiy_8,biy_8,ciy_8,
$ LDIST_DIM,Vspng_nk,trp_12emax,trp_22emax,G_nj)
endif
*
* TLM
* ---
call vspng_del2
(F_ip,
$ Opr_opsxp0_8(G_ni+1),Opr_opsyp0_8(G_nj+1),
$ aix_8,bix_8,cix_8,dix_8,aiy_8,biy_8,ciy_8,
$ LDIST_DIM,Vspng_nk,trp_12emax,trp_22emax,G_nj)
*
* TRAJECTORY
* ----------
if (V4dr_redotr_L) then
call vspng_del2
(F_itm,
$ Opr_opsxp0_8(G_ni+1),Opr_opsyp0_8(G_nj+1),
$ aix_8,bix_8,cix_8,dix_8,aiy_8,biy_8,ciy_8,
$ LDIST_DIM,Vspng_nk,trp_12emax,trp_22emax,G_nj)
endif
*
* TLM
* ---
call vspng_del2
(F_it,
$ Opr_opsxp0_8(G_ni+1),Opr_opsyp0_8(G_nj+1),
$ aix_8,bix_8,cix_8,dix_8,aiy_8,biy_8,ciy_8,
$ LDIST_DIM,Vspng_nk,trp_12emax,trp_22emax,G_nj)
*
!$omp parallel
*
if (V4dr_redotr_L) then
!$omp do
do k = 2, Vspng_nk
do j = 1, l_nj
do i = 1, l_ni
*
* TRAJECTORY
* ----------
F_itm(i,j,k) = F_ipm(i,j,k) + Cstvr_fistr_8(k) + F_xfis(i,j)
*
end do
end do
end do
!$omp enddo
endif
*
!$omp do
do k = 2, Vspng_nk
do j = 1, l_nj
do i = 1, l_ni
*
* TLM
* ---
F_it(i,j,k) = F_ip(i,j,k)
*
end do
end do
end do
!$omp enddo
*
* updating phi' at the top
*
if (V4dr_redotr_L) then
!$omp do
do j = 1, l_nj
do i = 1, l_ni
*
* TRAJECTORY
* ----------
F_ipm(i,j,1) = F_itm(i,j,1) - Cstvr_fistr_8(1) - F_xfis(i,j)
*
end do
end do
!$omp enddo
endif
*
!$omp do
do j = 1, l_nj
do i = 1, l_ni
*
* TLM
* ---
F_ip(i,j,1) = F_it(i,j,1)
*
end do
end do
!$omp enddo
*
* Save pi' for the diffusion of q'
* ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
if ( (.not.Schm_hydro_L) .and. Schm_difqp_L ) then
!$omp do
do k = 1, Vspng_nk
do j = 1, l_nj
do i = 1, l_ni
*
* TRAJECTORY
* ----------
work_m(i,j,k) = F_pipm(i,j,k)
*
* TLM
* ---
work(i,j,k) = F_pip(i,j,k)
*
end do
end do
end do
!$omp enddo
endif
!$omp end parallel
*
* The hydrostatic pressure: pi'
* ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
*
* TRAJECTORY
* ----------
call vspng_del2
(F_pipm,
$ Opr_opsxp0_8(G_ni+1),Opr_opsyp0_8(G_nj+1),
$ aix_8,bix_8,cix_8,dix_8,aiy_8,biy_8,ciy_8,
$ LDIST_DIM,Vspng_nk,trp_12emax,trp_22emax,G_nj)
*
* TLM
* ---
call vspng_del2
(F_pip,
$ Opr_opsxp0_8(G_ni+1),Opr_opsyp0_8(G_nj+1),
$ aix_8,bix_8,cix_8,dix_8,aiy_8,biy_8,ciy_8,
$ LDIST_DIM,Vspng_nk,trp_12emax,trp_22emax,G_nj)
*
*****************************
* 3. Nonyhydrostatic model *
*****************************
*
if ( .not. Schm_hydro_L ) then
*
* Vertical wind (physical)
* ~~~~~~~~~~~~~~~~~~~~~~~~
*
* TRAJECTORY
* ----------
if (V4dr_redotr_L) then
call vspng_del2
(F_wm,
$ Opr_opsxp0_8(G_ni+1),Opr_opsyp0_8(G_nj+1),
$ aix_8,bix_8,cix_8,dix_8,aiy_8,biy_8,ciy_8,
$ LDIST_DIM,Vspng_nk,trp_12emax,trp_22emax,G_nj)
endif
*
* TLM
* ---
call vspng_del2
(F_w,
$ Opr_opsxp0_8(G_ni+1),Opr_opsyp0_8(G_nj+1),
$ aix_8,bix_8,cix_8,dix_8,aiy_8,biy_8,ciy_8,
$ LDIST_DIM,Vspng_nk,trp_12emax,trp_22emax,G_nj)
*
if ( Schm_difqp_L ) then ! q' & related variables
*
* TRAJECTORY
* ----------
if (V4dr_redotr_L) then
call vspng_del2
(F_mum,
$ Opr_opsxp0_8(G_ni+1),Opr_opsyp0_8(G_nj+1),
$ aix_8,bix_8,cix_8,dix_8,aiy_8,biy_8,ciy_8,
$ LDIST_DIM,Vspng_nk,trp_12emax,trp_22emax,G_nj)
endif
*
* TLM
* ---
call vspng_del2
(F_mu,
$ Opr_opsxp0_8(G_ni+1),Opr_opsyp0_8(G_nj+1),
$ aix_8,bix_8,cix_8,dix_8,aiy_8,biy_8,ciy_8,
$ LDIST_DIM,Vspng_nk,trp_12emax,trp_22emax,G_nj)
*
* TRAJECTORY
* ----------
if (V4dr_redotr_L) then
call vspng_del2
(F_mulm,
$ Opr_opsxp0_8(G_ni+1),Opr_opsyp0_8(G_nj+1),
$ aix_8,bix_8,cix_8,dix_8,aiy_8,biy_8,ciy_8,
$ LDIST_DIM,Vspng_nk,trp_12emax,trp_22emax,G_nj)
endif
*
* TLM
* ---
call vspng_del2
(F_mul,
$ Opr_opsxp0_8(G_ni+1),Opr_opsyp0_8(G_nj+1),
$ aix_8,bix_8,cix_8,dix_8,aiy_8,biy_8,ciy_8,
$ LDIST_DIM,Vspng_nk,trp_12emax,trp_22emax,G_nj)
*
* Indirect diffusion of q'
* ~~~~~~~~~~~~~~~~~~~~~~~~
!$omp parallel do
do k=1, Vspng_nk
do j=1, l_nj
do i=1, l_ni
*
* TLM
* ---
work(i,j,k) = ( dble(work (i,j,k)) )
% *exp( dble(F_qpm(i,j,k)) ) +
% ( Geomg_z_8(k) + dble(work_m(i,j,k)) )
% *exp( dble(F_qpm(i,j,k)) ) * dble(F_qp(i,j,k))
*
* TRAJECTORY
* ----------
work_m(i,j,k) = ( Geomg_z_8(k) + dble(work_m(i,j,k)) )
% *exp( dble(F_qpm(i,j,k)) )
*
end do
end do
end do
!$omp end parallel do
*
* TRAJECTORY
* ----------
call vspng_del2
(work_m,
$ Opr_opsxp0_8(G_ni+1),Opr_opsyp0_8(G_nj+1),
$ aix_8,bix_8,cix_8,dix_8,aiy_8,biy_8,ciy_8,
$ LDIST_DIM,Vspng_nk,trp_12emax,trp_22emax,G_nj)
* TLM
* ---
call vspng_del2
(work,
$ Opr_opsxp0_8(G_ni+1),Opr_opsyp0_8(G_nj+1),
$ aix_8,bix_8,cix_8,dix_8,aiy_8,biy_8,ciy_8,
$ LDIST_DIM,Vspng_nk,trp_12emax,trp_22emax,G_nj)
*
!$omp parallel do
do k=1, Vspng_nk
do j=1, l_nj
do i=1, l_ni
*
* TRAJECTORY
* ----------
F_qpm(i,j,k) = log( dble(work_m(i,j,k))
% /( Geomg_z_8(k) + dble(F_pipm(i,j,k)) ) )
*
* TLM
* ---
F_qp(i,j,k) = (dble(work (i,j,k))*( Geomg_z_8(k) + dble(F_pipm(i,j,k))) -
% dble(work_m(i,j,k))*( dble(F_pip (i,j,k))) )
% /(dble(work_m(i,j,k))*( Geomg_z_8(k) + dble(F_pipm(i,j,k))) )
*
end do
end do
end do
!$omp end parallel do
*
endif
endif
*
if (Vspng_uvwdt_L) Vspng_nk = nkspng
*
endif
*
* __________________________________________________________________
*
return
end