!-------------------------------------- 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 simpsas(na_indic,na_dim,da_u,da_f,da_gradF) 1,19
use oda_shared
, only : dg_ubar
use mod4dv
, only : l4dvar
use modstag
, only : lstagwinds
implicit none
! Argument declarations
integer :: na_dim ! Dimension of the control vector in observation error covariances space
integer :: na_indic ! =1 No action taken; =4 Both F(u) and its gradient are computed.
real*8 :: da_F ! Objective function of the PSAS algorithm
real*8, dimension(na_dim) :: da_gradF ! Gradient of the PSAS objective funtion
real*8, dimension(na_dim) :: da_u ! Control variable in observation error covariances space
!
! Purpose: Implement the PSAS solver as described in
! Courtier, 1997, Dual formulation of four-dimentional variational assimilation,
! Q.J.R., pp2449-2461.
!
! Author : Simon Pellerin *ARMA/MSC January 2005
!
! Revisions:
!
#include "comlun.cdk"
#include "comdim.cdk"
#include "comcva.cdk"
#include "comoabdy.cdk"
! Local declarations
integer :: nl_ilev, nl_err
real*8, dimension(na_dim) :: dl_u,dl_d
if (na_indic .eq. 1 .or. na_indic .eq. 4) call tmg_stop (21)
call tmg_start(31,'SIMPSAS')
if (na_indic .ne. 1) then ! No action taken if na_indic == 1
nsim3d = nsim3d + 1
write(nulout,*) 'Entering simpsas for simulation ',nsim3d
dl_u = da_u + dg_ubar !dg_ubar = sum(u) of previous outer-loops
! Computation of uTu/2 term
da_F = dot_product(dl_u,dl_u)/2.d0
call tmg_start(2,'OBS-OPER') !
call oda_u2cma
(ncmoma,da_u,na_dim) ! transfer u in ncmoma
call oda_sqrtRm1
(ncmomi,ncmoma) ! sqrt(R-1)u
CALL TRANSFER
('ZOB0')
call oda_HT
! OMI -> Gomobs
call tmg_stop(2)
CALL TRANSFER
('ZGD0')
call oda_LT
! Gomobs -> GD
CALL TRANSFER
('ZSP0')
vazx = 0.d0
call oda_sqrtBT
(vazx,nvadim) ! GD -> vazx
call oda_sqrtB
(vazx,nvadim) ! vazx -> GD
call oda_L
! GD -> Gomobs
! Computation 'R -1/2 HBHT R-1/2 u' and store it for F and gradF computations
!
call tmg_start(2,'OBS-OPER') !
! Compute R-1/2 d and store it in ncmomi (d is store in ncmvar)
!
call oda_sqrtRm1
(ncmomi,ncmvar)
call oda_cma2u
(ncmomi,dl_d,na_dim) ! transfer ncmomi to dl_d (normalized inovation)
call oda_H
! Modify NCMOMA ; HBHT R-1/2 u ; Gomobs -> OMA
call oda_res
! HBHT R-1/2 u - d = Hdx - d -> OMA
! Computation of the normalized residual : sqrt(R-1) [Hdx -d]
call oda_sqrtRm1
(ncmoma,ncmoma) ! Modify NCMOMA : R-1/2 [HBHT R-1/2 u - d]
! da_gradF = oda_H(tg_y) / tg_y%rstd
call oda_cma2u
(ncmoma,da_gradF,na_dim) ! oma -> da_gradF
call tmg_stop(2) !
! at this point da_gradF contain R -1/2 [HBHT R-1/2 u - d]
! last term of objective function : F(u) = 1/2uTu + 1/2 uT {R-1/2[HBHT R-1/2 u - d] - R-1/2d]
da_F = da_F &
+ dot_product(dl_u,da_gradF-dl_d)/2.d0
! Gradient computation : gradF = u + R-1/2 (HBHT R-1/2 u -d)
da_gradF = dl_u + da_gradF
endif
call tmg_stop(31)
if (na_indic .eq. 1 .or. na_indic .eq. 4) call tmg_start(21,'QN')
end subroutine simpsas