!-------------------------------------- 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_testtlm - Verify the tangent linear model / direct model
*
#include "model_macros_f.h"
*
subroutine v4d_testtlm 1,17
*
implicit none
*
*author
* M.Tanguay
*
*revision
* v2_10 - Tanguay M. - initial MPI version
* v3_02 - Mahfouf J.-F. - add V4dg_anincr_L
* v3_03 - Tanguay M. - Adjoint Lam and NoHyd configuration
* v3_11 - Tanguay M. - Introduce V4dg_oktrcv_L
* - Add diagnostics of HU
* v3_31 - Tanguay M. - Control BC
*
*object
* see id section
*
*arguments
* none
*
*implicits
#include "glb_ld.cdk"
#include "lun.cdk"
#include "v4dg.cdk"
#include "v4dc.cdk"
#include "lctl.cdk"
#include "schm.cdk"
#include "ptopo.cdk"
#include "vt1.cdk"
#include "ind.cdk"
#include "tr3d.cdk"
#include "v4dg_bc.cdk"
*
*
integer vmmlod,vmmget,vmmuld
external vmmlod,vmmget,vmmuld
*
integer n,nn,lambda,lima,limb,wrt_out,indic,i0,j0,in,jn,key1(2),err,
% kk,istart,status
real pj
*
* -----------------------------------------------------------------
* Parameters for list of stations and variables where TLM is tested
* -----------------------------------------------------------------
integer NSTA_H,NVAR_H,NSTA_NH,NVAR_NH,NSTA_TOT,NVAR_TOT,NSTA_TR,NVAR_TR
*
* --------------------------------------
* Range of exponents used in TLM testing
* --------------------------------------
integer RANGMAX
*
real, dimension (:,:,:), allocatable:: tab
integer, dimension (:) , allocatable:: ni,nj,nk,nprev,npos,ix,jx,kx,lx
* _________________________________________________________________________
*
if (Lun_out.gt.0) then
write(Lun_out,fmt=
$ '('' TEST OF TANGENT LINEAR with CONF = '',I6)')V4dg_conf
write(Lun_out,fmt=
$ '('' --------------------------------------'')')
endif
*
* Initializations
* ---------------
RANGMAX = 5
*
NSTA_H = 10
NVAR_H = 4
*
if(.not.Schm_hydro_L) then
NSTA_NH = 2
NVAR_NH = 1
else
NSTA_NH = 0
NVAR_NH = 0
endif
*
* Tracers
* -------
NSTA_TR = 3
NVAR_TR = Tr3d_ntr
*
NSTA_TOT= NSTA_H + NSTA_NH + NSTA_TR
NVAR_TOT= NVAR_H + NVAR_NH + NVAR_TR
*
* Table Allocations
* -----------------
allocate (tab(8,RANGMAX+1,NSTA_TOT),
% ni(NVAR_TOT),nj(NVAR_TOT),nk(NVAR_TOT),nprev(NVAR_TOT),
% npos(NSTA_TOT),ix(NSTA_TOT),jx(NSTA_TOT),kx(NSTA_TOT),lx(NSTA_TOT))
*
* Select array without pilot region when G_lam
* NOTE: Should be valid for all variables
* --------------------------------------------
i0 = 1+pil_w
in = l_niu-pil_e
j0 = 1+pil_s
jn = l_njv-pil_n
*
* ------------------------------------------------------
* Set list of stations and variables where TLM is tested
* ------------------------------------------------------
*
* ----------------------------------------------------
* Predictives variables available for hydrostatic runs
* (1=UU,2=VV,3=4T,4=4S)
* ----------------------------------------------------
*
* Dimensions
* ----------
ni(1) = l_niu
ni(2) = l_ni
ni(3) = l_ni
ni(4) = l_ni
*
nj(1) = l_nj
nj(2) = l_njv
nj(3) = l_nj
nj(4) = l_nj
*
nk(1) = l_nk
nk(2) = l_nk
nk(3) = l_nk
nk(4) = 1
*
* Set positions
* -------------
do nn = 1, NSTA_H
*
* -------------------------------------
* LX= Range of a predictive variable
* IX= I indice of a predictive variable
* JX= J indice of a predictive variable
* KX= K indice of a predictive variable
* -------------------------------------
lx(nn) = mod(nn,NVAR_H) + 1
ix(nn) = (float(nn-1)/NSTA_H)*(in-i0+1) + i0
jx(nn) = (float(nn-1)/NSTA_H)*(jn-j0+1) + j0
kx(nn) = (float(nn-1)/NSTA_H)*l_nk + 1
*
* LX=4 (4S) is Surface Variable
* -----------------------------
if(lx(nn).eq.4) kx(nn)=1
*
end do
*
* --------------------------------------------------------
* Predictives variables available for non hydrostatic runs
* (5=FIP)
* --------------------------------------------------------
*
if(.not.Schm_hydro_L) then
*
* Dimensions
* ----------
ni(5) = l_ni
nj(5) = l_nj
nk(5) = l_nk
*
* Set positions
* -------------
do nn = NSTA_H+1, NSTA_H + NSTA_NH
*
* -------------------------------------
* LX= Range of a predictive variable
* IX= I indice of a predictive variable
* JX= J indice of a predictive variable
* KX= K indice of a predictive variable
* -------------------------------------
lx(nn) = mod(nn,NVAR_NH) + 1
lx(nn) = lx(nn) + NVAR_H
ix(nn) = (float(nn-(NSTA_H+1))/NSTA_NH)*(in-i0+1) + i0
jx(nn) = (float(nn-(NSTA_H+1))/NSTA_NH)*(jn-j0+1) + j0
kx(nn) = (float(nn-(NSTA_H+1))/NSTA_NH)*l_nk + 1
*
end do
*
endif
*
* -------------------------------
* Predictives 3D tracer variables
* -------------------------------
*
istart = NVAR_H + NVAR_NH
*
* Dimensions
* ----------
do kk = 1,NVAR_TR
ni(kk+istart) = l_ni
nj(kk+istart) = l_nj
nk(kk+istart) = l_nk
end do
*
* Set positions
* -------------
do nn = NSTA_H+NSTA_NH+1, NSTA_TOT
*
* -------------------------------------
* LX= Range of a predictive variable
* IX= I indice of a predictive variable
* JX= J indice of a predictive variable
* KX= K indice of a predictive variable
* -------------------------------------
*
C Get HU only
C -----------
C lx(nn) = mod(nn,NVAR_TR) + 1
lx(nn) = 1
*
lx(nn) = lx(nn) + NVAR_H + NVAR_NH
ix(nn) = (float(nn-(NSTA_H+NSTA_NH+1))/NSTA_TR)*(in-i0+1) + i0
jx(nn) = (float(nn-(NSTA_H+NSTA_NH+1))/NSTA_TR)*(jn-j0+1) + j0
kx(nn) = (float(nn-(NSTA_H+NSTA_NH+1))/NSTA_TR)*l_nk + 1
*
end do
*
* Block Position of previous predictive variable
* ----------------------------------------------
nprev(1) = 0
do nn = 2,NVAR_TOT
nprev(nn) = nprev(nn-1) + ni(nn-1) * nj(nn-1) * nk(nn-1)
end do
*
* --------------------------------------
* Read given analysis of model variables
* --------------------------------------
if( Lun_out.gt.0 ) then
write(Lun_out, fmt='(//''-------------------'')')
write(Lun_out, fmt='( ''READ GIVEN ANALYSIS'')')
write(Lun_out, fmt='( ''-------------------'')')
endif
*
V4dg_part = 2
call indata
()
*
* ------------------------------------------------------------
* Initialize starting control var. for REFERENCE (NLM) and TLM
* ------------------------------------------------------------
*
* Initialize REFERENCE (NLM) initial control var. in V4dc_ycv
* from model var.
* -----------------------------------------------------------
call v4d_cainin
(V4dc_ncv,V4dc_ycv)
*
* Use a hydrostatic phi' when .not.Schm_hydro_L
* ---------------------------------------------
if(.not.Schm_hydro_L) call v4d_phydro
(V4dc_ncv,V4dc_ycv)
*
if ( V4dg_anincr_L ) then
*
* Read perturbed (NLM) initial conditions for trial field
* to create realistic perturbations in V4dc_xcv
* -------------------------------------------------------
if( Lun_out.gt.0 ) then
write(Lun_out, fmt='(//''----------------'')')
write(Lun_out, fmt='( ''READ GIVEN TRIAL'')')
write(Lun_out, fmt='( ''----------------'')')
endif
*
call v4d_rdtrial
()
*
* Get fields in memory
* --------------------
key1(1) = VMM_KEY(ut1)
key1(2) = VMM_KEY(vt1)
err = vmmlod(key1,2)
err = VMM_GET_VAR(ut1)
err = VMM_GET_VAR(vt1)
*
* Associate with Ind
* ------------------
Ind_u_ = ut1_
Ind_v_ = vt1_
*
* Convert wind images to true winds
* ---------------------------------
call v4d_img2uv
()
*
err = vmmuld(-1,0)
*
* Initialize V4dc_xcv from trial field
* ------------------------------------
call v4d_cainin
(V4dc_ncv,V4dc_xcv)
*
* Initialize TLM initial control var. in V4dc_xcv
* -----------------------------------------------
do n = 1,V4dc_ncv
V4dc_xcv(n) = V4dc_ycv(n) - V4dc_xcv(n)
end do
*
else
*
* Keep NLM initial control var. in V4dc_xcv
* NOTE: Used later to define TLM initial control var.
* ---------------------------------------------------
do n = 1,V4dc_ncv
V4dc_xcv(n) = V4dc_ycv(n)
end do
*
endif
*
if( V4dg_conf.ne.510 ) then
* -------------------------------------------
* Keep NLM initial control var. in V4dc_wkmin
* -------------------------------------------
*
do n = 1,V4dc_ncv
V4dc_wkmin( 4*V4dc_ncv + n ) = V4dc_ycv(n)
end do
*
else
* -------------------------------------------
* Keep NLM initial control var. in V4dc_wkmin
* -------------------------------------------
*
do n = 1,V4dc_ncv
V4dc_wkmin(n)= V4dc_ycv(n)
end do
*
endif
*
* -------------------------------
* Run REFERENCE (NLM) integration
* -------------------------------
if( Lun_out.gt.0 .and. Ptopo_myproc.eq.0 ) then
write(Lun_out, fmt='(//''---------------------'')')
write(Lun_out, fmt='( ''REFERENCE INTEGRATION'')')
write(Lun_out, fmt='( ''---------------------'')')
endif
*
* Set status of the integration (REFERENCE integration)
* -----------------------------------------------------
V4dg_status = 5
*
* Set over-riding switch for dynout and blocstat
* ----------------------------------------------
V4dg_output_L = .true.
*
* Call 4D-Var simulator with NLM integration without gradient
* -----------------------------------------------------------
indic = 99
call v4d_simul
(indic,V4dc_ncv,V4dc_ycv,pj,V4dc_gcv)
*
* Initialize V4dc_ycv (Image winds) from NLM final model var.
* -----------------------------------------------------------
call v4d_cainin
(V4dc_ncv,V4dc_ycv)
*
* -----------------------------------------
* Keep NLM final control var. in V4dc_wkmin
* -----------------------------------------
if( V4dg_conf.ne.510 ) then
*
do n = 1,V4dc_ncv
V4dc_wkmin(n) = V4dc_ycv(n)
end do
*
endif
*
* Get fields in memory
* --------------------
key1(1) = VMM_KEY(ut1)
key1(2) = VMM_KEY(vt1)
err = vmmlod(key1,2)
err = VMM_GET_VAR(ut1)
err = VMM_GET_VAR(vt1)
*
* Associate with Ind
* ------------------
Ind_u_ = ut1_
Ind_v_ = vt1_
*
* Convert wind images to true winds
* ---------------------------------
call v4d_img2uv
()
*
err = vmmuld(-1,0)
*
* Initialize V4dc_ycv (True winds)) from NLM final model var.
* -----------------------------------------------------------
call v4d_cainin
(V4dc_ncv,V4dc_ycv)
*
* ------------------------------
* Run TANGENT LINEAR integration
* ------------------------------
if( Lun_out.gt.0 .and. Ptopo_myproc.eq.0 ) then
write(Lun_out, fmt='(//''--------------------------'')')
write(Lun_out, fmt='( ''TANGENT LINEAR INTEGRATION'')')
write(Lun_out, fmt='( ''--------------------------'')')
endif
*
* Define TLM control var. with controlled size perturbations
* ----------------------------------------------------------
if ( .not.V4dg_anincr_L ) then
*
do n = 1,V4dc_ncv
V4dc_xcv(n) = V4dc_ycv(n) - V4dc_xcv(n)
end do
*
endif
*
* Use a TLM hydrostatic phi' when .not.Schm_hydro_L
* -------------------------------------------------
if(.not.Schm_hydro_L) call v4d_phydro_tl
(V4dc_ncv,V4dc_xcv)
*
* Zero pilot region of perturbation fields when G_lam
* ---------------------------------------------------
if(G_lam.and.V4dg_bc_variant.eq.0) call v4d_zeropilot
(V4dc_ncv,V4dc_xcv)
*
* --------------------------------------------
* Keep TLM initial control var. in V4dc_wkmin
* --------------------------------------------
do n = 1,V4dc_ncv
V4dc_wkmin(V4dc_ncv+n) = V4dc_xcv(n)
end do
*
* Set status of the integration (TLM integration)
* -----------------------------------------------
V4dg_status = 10
*
* Set over-riding switch for dynout and blocstat
* ----------------------------------------------
V4dg_output_L = .true.
*
* Call 4D-Var simulator with TLM model without gradient
* -----------------------------------------------------
V4dg_tlm_L = .true.
indic = 99
call v4d_simul
(indic,V4dc_ncv,V4dc_xcv,pj,V4dc_gcv)
*
* Reset 4D-Var simulator with NLM model
* -------------------------------------
V4dg_tlm_L = .false.
*
* Initialize V4dc_ycv from TLM final model variables
* --------------------------------------------------
call v4d_cainin
(V4dc_ncv,V4dc_ycv)
*
* -----------------------------------------
* Keep TLM final control var. in V4dc_wkmin
* -----------------------------------------
if( V4dg_conf.ne.510 ) then
*
do n = 1,V4dc_ncv
V4dc_wkmin( 2*V4dc_ncv+n ) = V4dc_ycv(n)
end do
*
endif
*
* Deactivate writing on TRAJECTORY VMM
* ------------------------------------
V4dg_oktr_L = .false.
*
* Deactivate writing on TRAJECTORY Conversion
* -------------------------------------------
V4dg_oktrcv_L = .false.
*
if( V4dg_conf.eq.510 ) then
* --------------------------------------------
* Do only CURRENT (NLM) integration and RETURN
* --------------------------------------------
*
if( Lun_out.gt.0 .and. Ptopo_myproc.eq.0 ) then
write(Lun_out, fmt='(//''-------------------'')')
write(Lun_out, fmt='( ''CURRENT INTEGRATION'')')
write(Lun_out, fmt='( ''-------------------'')')
endif
*
* Set CURRENT initial control var. as REFERENCE + TLM
* ---------------------------------------------------
do n = 1,V4dc_ncv
V4dc_ycv(n) = V4dc_wkmin(n) + V4dc_wkmin(V4dc_ncv+n)
end do
*
* Set status of the integration (CURRENT integration)
* ---------------------------------------------------
V4dg_status = 0
*
* Set over-riding switch for dynout and blocstat
* ----------------------------------------------
V4dg_output_L = .true.
*
* Call 4D-Var simulator with NLM integration without gradient
* -----------------------------------------------------------
indic = 99
call v4d_simul
(indic,V4dc_ncv,V4dc_ycv,pj,V4dc_gcv)
*
else
* ------------------------------------------------
* Ratio TLM test as in Thepaut and Courtier (1991)
* ------------------------------------------------
*
if( Lun_out.gt.0 .and. Ptopo_myproc.eq.0 ) then
write(Lun_out, fmt='(/''-------------------------------------'')')
write(Lun_out, fmt='( ''Ratio TLM test (THEPAUT AND COURTIER)'')')
write(Lun_out, fmt='( ''-------------------------------------''/)')
endif
*
lima = 0
limb = RANGMAX
*
* Loop over lambda ( pert. = 10**(-lambda) * TLM pert. )
* ------------------------------------------------------
do lambda = limb,lima,-1
*
if( Lun_out.gt.0 ) then
write(Lun_out, fmt='(/''--------------------------------------'')')
write(Lun_out, fmt='( ''CURRENT INTEGRATION FOR LAMBDA ='',I6)') lambda
write(Lun_out, fmt='( ''--------------------------------------''/)')
endif
*
* Set CURRENT initial control var. as REFERENCE + 10**(-lambda)*TLM
* -----------------------------------------------------------------
do n = 1,V4dc_ncv
V4dc_ycv(n) = V4dc_wkmin( 4*V4dc_ncv+n) +
% V4dc_wkmin( V4dc_ncv+n) * 10.**(-lambda)
end do
*
* -----------------------------------------------
* Keep CURRENT initial control var. in V4dc_wkmin
* -----------------------------------------------
do n = 1,V4dc_ncv
V4dc_wkmin(3*V4dc_ncv+n) = V4dc_ycv(n)
end do
*
* Set status of the integration (CURRENT integration)
* ---------------------------------------------------
V4dg_status = 0
*
* Set over-riding switch for dynout and blocstat
* ----------------------------------------------
V4dg_output_L = .false.
if(lambda.eq.lima) V4dg_output_L = .true.
*
* Cancel write(Lun_out)
* ---------------------
wrt_out = Lun_out
if(lambda.ne.lima) Lun_out = -99
*
* Call 4D-Var simulator with NLM integration without gradient
* -----------------------------------------------------------
indic = 99
call v4d_simul
(indic,V4dc_ncv,V4dc_ycv,pj,V4dc_gcv)
*
* Reset write(Lun_out)
* --------------------
Lun_out = wrt_out
*
* Initialize V4dc_ycv from NLM final model var.
* ---------------------------------------------
call v4d_cainin
(V4dc_ncv,V4dc_ycv)
*
* ----------------------------------------------
* Keep CURRENT final control var. in V4dc_wkmin
* ----------------------------------------------
do n = 1,V4dc_ncv
V4dc_wkmin(5*V4dc_ncv + n) = V4dc_ycv(n)
end do
*
* -------------------------------------------------------------
* Computation of CURRENT (NLM) - REFERENCE (NLM) divided by TLM
* at final state for prescribed locations
* -------------------------------------------------------------
do nn = 1, NSTA_TOT
*
* Set location according to station and variable
* ----------------------------------------------
npos(nn)= nprev(lx(nn)) + (kx(nn)-1) * nj(lx(nn)) * ni(lx(nn)) +
% (jx(nn)-1) * ni(lx(nn)) + ix(nn)
*
n = npos(nn)
*
* Do computations
* ---------------
V4dc_wkmin(6*V4dc_ncv+n) = V4dc_wkmin(5*V4dc_ncv+n)-V4dc_wkmin(n)
*
if(V4dc_wkmin(2*V4dc_ncv+n).ne.0.) then
V4dc_wkmin(7*V4dc_ncv+n) = V4dc_wkmin(6*V4dc_ncv+n)/
% (V4dc_wkmin(2*V4dc_ncv+n)*10.**(-lambda))
else
V4dc_wkmin(7*V4dc_ncv+n) =-5555555
endif
*
* Store results
* -------------
tab(1,lambda+1,nn) = V4dc_wkmin(7*V4dc_ncv+n)
tab(2,lambda+1,nn) = V4dc_wkmin(6*V4dc_ncv+n)
tab(3,lambda+1,nn) = V4dc_wkmin(2*V4dc_ncv+n)
tab(4,lambda+1,nn) = V4dc_wkmin(5*V4dc_ncv+n)
tab(5,lambda+1,nn) = V4dc_wkmin(3*V4dc_ncv+n)
tab(6,lambda+1,nn) = V4dc_wkmin(0*V4dc_ncv+n)
tab(7,lambda+1,nn) = V4dc_wkmin(4*V4dc_ncv+n)
tab(8,lambda+1,nn) = V4dc_wkmin(1*V4dc_ncv+n)
*
* ENDDO NN
* --------
end do
*
* ENDDO LAMBDA
* ------------
end do
*
do nn = 1, NSTA_TOT
*
* Recall location according to station and variable
* -------------------------------------------------
n = npos(nn)
*
* Print results
* -------------
if( Lun_out.gt.0 ) then
*
write(Lun_out,3000) n,Ptopo_myproc
if(lx(nn).eq.1) write(Lun_out,3001) ix(nn),jx(nn),kx(nn)
if(lx(nn).eq.2) write(Lun_out,3002) ix(nn),jx(nn),kx(nn)
if(lx(nn).eq.3) write(Lun_out,3003) ix(nn),jx(nn),kx(nn)
if(lx(nn).eq.4) write(Lun_out,3004) ix(nn),jx(nn)
if(lx(nn).eq.5.and..not.Schm_hydro_L) write(Lun_out,3005) ix(nn),jx(nn),kx(nn)
if(lx(nn).eq.5.and. Schm_hydro_L) write(Lun_out,3006) ix(nn),jx(nn),kx(nn)
if(lx(nn).eq.6.and..not.Schm_hydro_L) write(Lun_out,3006) ix(nn),jx(nn),kx(nn)
*
do lambda = lima,limb
write(Lun_out,4000) -lambda,(tab(indic,lambda+1,nn),indic=1,8)
end do
*
endif
*
* ENDDO NN
* --------
end do
*
endif
*
deallocate( tab, STAT=status )
*
3000 format(
+//,'TEST OF THE TANGENT LINEAR AT POINT = ',I6,' PROCESS = ',I4,
+ /,'==========================================================')
3001 format(
+'U WIND AT (I,J,K) = (',I5,',',I5,',',I5,')',
+ /,'==========================================================',
+//)
3002 format(
+'V WIND AT (I,J,K) = (',I5,',',I5,',',I5,')',
+ /,'==========================================================',
+//)
3003 format(
+'T PRIME AT (I,J,K) = (',I5,',',I5,',',I5,')'
+ /,'==========================================================',
+//)
3004 format(
+'S AT (I,J) = (',I5,',',I5,')',
+ /,'==========================================================',
+//)
3005 format(
+'FP PRIME AT (I,J,K) = (',I5,',',I5,',',I5,')'
+ /,'==========================================================',
+//)
3006 format(
+'HU PRIME AT (I,J,K) = (',I5,',',I5,',',I5,')'
+ /,'==========================================================',
+//)
*
4000 format(
+'LAMBDA =',I3,' RATIO = ',E14.8,2X, 7(e11.5,1x))
*
return
end