!-------------------------------------- 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_testgrd - Routine which tests the correctness of the gradient
*
#include "model_macros_f.h"
*

      subroutine v4d_testgrd (Ndim,F_px,F_pgraj,F_py,F_pj0,F_start,F_range) 2,4
*
      implicit none
*
      integer Ndim,F_range
      real F_px(Ndim),F_pgraj(Ndim),F_py(Ndim),F_pj0,F_start
*
*author
*     M.Tanguay
*
*revision
* v2_10 - Tanguay M.        - initial MPI version
* v3_11 - Tanguay M.        - Introduce V4dg_oktrcv_L 
*
*object
*     Evaluate the correctness of the gradient
*     using the Taylor formula
*
*arguments
*  Name        I/O                 Description
*----------------------------------------------------------------
* Ndim         I                   Dimension 
* F_px         I                   Control variable at initial time
* F_pgraj      O                   Gradient at initial time
* F_py         -                   Work control variable 
* F_pj0        O                   Cost function value
* F_start      I                   Starting factor 
* F_range      I                   Range of exponents                  
*-----------------------------------------------------------------
*
*implicits
#include "glb_ld.cdk"
#include "lun.cdk"
#include "v4dg.cdk"
*
      integer k,indic,jt,wrt_out
      logical out_L
*
      real xend,alp,pjf,test
      real*8 gnorm_8
*
*     ---------------------------------------------------
*     Range of exponents in alp = 10**(-exponent)*F_start 
*     ---------------------------------------------------
      integer RANGMAX
      parameter ( RANGMAX = 10 ) 
      real tab(RANGMAX+1,3)
*     ______________________________________________________
*
      if(F_range.gt.RANGMAX) call gefstop('v4d_testgrd')
*     ______________________________________________________
*
*     -----------------------------------------------------------------
*     Call 4D-Var simulator at F_px to evaluate the gradient in F_pgraj
*     -----------------------------------------------------------------
      indic = 4
      call v4d_simul (indic,Ndim,F_px,F_pj0,F_pgraj)
*
*     --------------------------------------
*     Computation of norm**2 of the gradient
*     --------------------------------------
      call v4d_scalpro (Ndim,F_pgraj,F_pgraj,gnorm_8)
*
*     ----------------
*     Perform the test
*     ----------------
*
      do jt=0,F_range
*
         alp = 10.**(-jt)*F_start
*
         do k = 1,Ndim
            F_py(k) = F_px(k) - alp*F_pgraj(k)
         end do
*
*        Set differently over-riding switch for dynout and blocstat 
*        ----------------------------------------------------------
         out_L         = V4dg_output_L 
         V4dg_output_L = .false.
*
*        Cancel write(Lun_out)
*        ---------------------
         wrt_out = Lun_out
         Lun_out = -99
*
*        Deactivate WRITE option on TRAJECTORY VMM WA file
*        -------------------------------------------------
         V4dg_oktr_L = .false.
*
*        Deactivate WRITE option on TRAJECTORY Conversion WA file
*        --------------------------------------------------------
         V4dg_oktrcv_L = .false.
*
*        Call 4D-Var simulator without gradient
*        --------------------------------------
         indic = 99
         call v4d_simul (indic,Ndim,F_py,pjf,F_pgraj)
*
*        Reactivate WRITE option on TRAJECTORY VMM WA file
*        -------------------------------------------------
         V4dg_oktr_L = .true.
*
*        Reactivate WRITE option on TRAJECTORY Conversion WA file
*        --------------------------------------------------------
         V4dg_oktrcv_L = .true.
*
*        Reset write(Lun_out)
*        --------------------
         Lun_out = wrt_out 
*
*        Reset over-riding switch for dynout and blocstat 
*        ------------------------------------------------
         V4dg_output_L = out_L 
*
*        Computes test and store result
*        ------------------------------
         test = (pjf-F_pj0)/(-alp*gnorm_8)
*
         tab(jt+1,1) = alp
         tab(jt+1,2) = pjf
         tab(jt+1,3) = test
*
      end do
*    
*     -----------
*     Diagnostics
*     -----------
      if(Lun_out.gt.0) then
*
         xend = F_start*10.**(-F_range)
         write(Lun_out,fmt=9100) F_start,xend,F_pj0,sqrt(sngl(gnorm_8))
         write(Lun_out,fmt=9200)
*
         do jt=0,F_range
            write(Lun_out,fmt=9201) jt,tab(jt+1,1),tab(jt+1,2),tab(jt+1,3)
         end do
            write(Lun_out,fmt=9202)
*
      endif
*
 9100 format(/,4X,'V4D_TESTGRD- The gradient is being tested for '
     %     ,E14.6,' <= ALPHA <= ',E14.6,/,17X,'Initial value of J(X):'
     %     ,1x,E14.8,4x,'Norm of GRAD J(X): ',E14.8)
 9200 format(/,5X,'I',10X,'ALPHA',12X,'J(X)',12X,'TEST')
 9201 format(3X,I3,4X,E14.6,4X,E14.8,4X,G14.8)
 9202 format(//)
*
      return
      end