!-------------------------------------- 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 simvar(na_indic,na_dim,da_v,da_J,da_gradJ) 6,20
  use oda_shared, only : dg_vbar
  use mod4dv, only : l4dvar
  use modstag, only : lstagwinds
  USE procs_topo
  implicit none
  ! Argument declarations
  integer :: na_dim ! Dimension of the control vector in forecast error coraviances space
  ! Value of na_indic
  ! Note: 1 and 4 are reserved values for call back from m1qn3.
  !       For direct calls use other value than 1 and 4.
  ! =1 No action taken; =4 Both J(u) and its gradient are computed.
  ! =2 Same as 4 (compute J and gradJ) but do not interrupt timer of the
  !    minimizer.
  ! =3 Compute Jo and gradJo only.
  integer :: na_indic 
  real*8  :: da_J ! Cost function of the Variational algorithm
  real*8, dimension(na_dim) :: da_gradJ ! Gradient of the Variational Cost funtion
  real*8, dimension(na_dim) :: da_v ! Control variable in forecast error covariances space
  !
  ! Purpose: Implement the Variational solver as described in
  ! Courtier, 1997, Dual formulation of four-dimentional variational assimilation,
  ! Q.J.R., pp2449-2461.
  !
  ! Author : Simon Pellerin *ARMA/MSC October 2005
  !          (Based on previous versions of evaljo.ftn, evaljg.ftn and evaljgns.ftn).
  !
  ! Revisions:
  ! L. Fillion, ARMA/EC, 5 Jun 2009. Introduce 1 Obs experiment.
  ! L. Fillion, ARMA/EC, 11 May 2010. Limit writing diagnostics to processor 0.
  !

#include "comlun.cdk"
#include "comvarqc.cdk"
#include "comoabdy.cdk"
#include "comdim.cdk"
#include "comcva.cdk"
  ! Local declaration
  integer :: nl_ilev, nl_err
  real*8, dimension(na_dim) :: dl_v
  real*8 :: dl_Jb, dl_Jo

  if (na_indic .eq. 1 .or. na_indic .eq. 4) call tmg_stop (21)

  call tmg_start(31,'SIMVAR')
  if (na_indic .ne. 1) then ! No action taken if na_indic == 1
     nsim3d = nsim3d + 1

     IF(myid == 0) write(nulout,*) 'Entering simvar for simulation ',nsim3d

     dl_v = da_v + dg_vbar      !dg_vbar = sum(v) of previous outer-loops

     ! Computation of vTv/2 term
     dl_Jb = dot_product(dl_v,dl_v)/2.d0

     call oda_sqrtB(da_v,na_dim)
     call oda_L

     call tmg_start(2,'OBS-OPER') !

     call oda_H ! Modify NCMOMA ! Hdx

     call oda_res              ! Modify NCMOMA : Hdx-d
     
     call oda_sqrtRm1(ncmoma,ncmoma)  ! Modify NCMOMA : R**-1/2*(Hdx-d)
     
     call oda_Jo               ! Store J in NCMOMI : R**-1(Hdx-d)**2
     
     IF (LVARQC) THEN
        ! Store modify J in NCMOMI : -ln((gamma-exp(J))/(gamma+1))
        call oda_qcv
     endif

     dl_Jo = 0.d0
     call oda_sumJo(dl_Jo)
     da_J = dl_Jb + dl_Jo
     if (na_indic .eq. 3) then
        da_J = dl_Jo
        IF(myid == 0) write(nulout,FMT='(6X,"SIMVAR:  JO = ",G23.16,6X)') dl_Jo
     else
        da_J = dl_Jb + dl_Jo
        IF(myid == 0) write(nulout,FMT='(6X,"SIMVAR:  Jb = ",G23.16,6X,"JO = ",G23.16,6X,"Jt = ",G23.16)') dl_Jb,dl_Jo,da_J
     endif

     call oda_sqrtRm1(ncmomi,ncmoma)

     IF (LVARQC) THEN
        call oda_qcvad
     endif

     CALL TRANSFER('ZOB0')
     call oda_HT
     
     call tmg_stop(2)

     CALL TRANSFER('ZGD0')
     call oda_LT

     CALL TRANSFER('ZSP0')
     da_gradJ = 0.d0
     call oda_sqrtBT(da_gradJ,na_dim)

     if (na_indic .ne. 3) then
        da_gradJ = dl_v + da_gradJ
     endif

  endif
  call tmg_stop(31)
  if (na_indic .eq. 1 .or. na_indic .eq. 4) call tmg_start(21,'QN')
end subroutine simvar