!-------------------------------------- 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 adw_trajex_ad - ADJ of adw_trajex_tl 
*
#include "model_macros_f.h"
*

      subroutine adw_trajex_ad ( F_xto,   F_yto,   F_xcto,  F_ycto, 3
     %                           F_zcto,  F_xctm,  F_yctm,  F_zctm,
     %                                    F_xctm_m,F_yctm_m,F_zctm_m,i0,in,j0,jn)
*
      implicit none
*
      real   F_xto (*), F_yto (*), F_xcto(*), F_ycto(*),
     %       F_zcto(*), F_xctm(*), F_yctm(*), F_zctm(*)
*
      real   F_xctm_m(*),F_yctm_m(*),F_zctm_m(*) 
      integer i0,in,j0,jn
*
*author
*     monique tanguay
*
*revision
* v2_31 - Tanguay M.        - initial MPI version
* v3_00 - Tanguay M.        - adapt to restructured adw_main 
* v3_01 - Tanguay M.        - correction minmax_m = +-1
* v3_11 - Tanguay M.        - AIXport+Opti+OpenMP for TLM-ADJ
*
*language
*     fortran 77
*
*object
*     see id section
*
*ADJ of
*arguments
*______________________________________________________________________
*              |                                                 |     |
* NAME         | DESCRIPTION                                     | I/O |
*--------------|-------------------------------------------------|-----|
*              |                                                 |     |
* F_xto        | upstream x positions at origin                  |  o  |
* F_yto        | upstream y positions at origin                  |  o  |
* F_xcto       | upstream x cartesian positions at origin        |  o  |
* F_ycto       | upstream y cartesian positions at origin        |  o  |
* F_zcto       | upstream z cartesian positions at origin        |  o  |
* F_xctm       | upstream x cartesian positions at mid-traj.     |  i  |
* F_yctm       | upstream y cartesian positions at mid-traj.     |  i  |
* F_zctm       | upstream z cartesian positions at mid-traj.     |  i  |
*______________|_________________________________________________|_____|
*
*implicits
#include "glb_ld.cdk"
#include "adw.cdk"
#include "dcst.cdk"
*
************************************************************************
      integer i,j,k, n, ij, nij
*
      real*8 prx_8, pry_8, prz_8, prdot2_8, r2pi_8, TWO_8, ZERO_8 
*
      parameter (TWO_8  = 2.0)
      parameter (ZERO_8 = 0.0)
*
      real*8 prdot2m_8
*
      real   prxcto_m, prycto_m, przcto_m, minmax, minmax_m
*
************************************************************************
      nij  = l_ni*l_nj
*
      r2pi_8 = TWO_8 * Dcst_pi_8
*
************************************************************************
*
*     Zero adjoint variables
*     ----------------------
      minmax   = ZERO_8
      prdot2_8 = ZERO_8
*
!$omp parallel private(n,ij,prx_8,pry_8,prz_8,prdot2_8,prdot2m_8,
!$omp&      prxcto_m,prycto_m,przcto_m,minmax,minmax_m)
!$omp do
      do k=1,l_nk
      do j=j0,jn
      do i=i0,in
      n = (k-1)*nij+((j-1)*l_ni) + i
*
      ij = mod( n-1, nij ) + 1
*
*     START REBUILD TRAJECTORY
*     ------------------------
      pry_8 = dble(Adw_cy2d_8(ij))
      prx_8 = dble(Adw_cx2d_8(ij)) * pry_8
      pry_8 = dble(Adw_sx2d_8(ij)) * pry_8
      prz_8 = dble(Adw_sy2d_8(ij))
*
*     TRAJECTORY
*     ----------
      prdot2m_8= 2.0 * ( prx_8 * dble(F_xctm_m(n)) +
     %                   pry_8 * dble(F_yctm_m(n)) +
     %                   prz_8 * dble(F_zctm_m(n)) )
*
      prxcto_m = prdot2m_8 * dble(F_xctm_m(n)) - prx_8
      prycto_m = prdot2m_8 * dble(F_yctm_m(n)) - pry_8
      przcto_m = prdot2m_8 * dble(F_zctm_m(n)) - prz_8
*
*     The following min statement is expanded as two IF blocks:
*     minmax_m = max(-1.,min(1.,przcto_m))
*
      minmax_m = przcto_m
*
*     END REBUILD TRAJECTORY
*     ----------------------
*
*     ADJ
*     ---
      if (przcto_m.ge.1.) then
              F_yto(n) = 0.
              minmax   = 0.
      elseif (przcto_m.le.-1.) then
              F_yto(n) = 0.
              minmax   = 0.
      else
C             minmax   = F_yto(n)/sqrt( 1.0-minmax_m**2 ) + minmax
              minmax   = F_yto(n)/sqrt( 1.0-minmax_m**2 )
              F_yto(n) = ZERO_8
      endif
      F_zcto(n) = minmax + F_zcto(n)
C     minmax    = ZERO_8
*
      F_ycto(n) = (  F_xto(n)*prxcto_m  )
     %              /(prxcto_m*prxcto_m + prycto_m*prycto_m)
     %              + F_ycto(n)
      F_xcto(n) = (  - prycto_m*F_xto(n))
     %              /(prxcto_m*prxcto_m + prycto_m*prycto_m)
     %              + F_xcto(n)
      F_xto(n)  = ZERO_8
*
      F_zctm(n) =     prdot2m_8 * dble(F_zcto(n))  + F_zctm(n)
*
      F_yctm(n) =     prdot2m_8 * dble(F_ycto(n))  + F_yctm(n)
*
      F_xctm(n) =     prdot2m_8 * dble(F_xcto(n))  + F_xctm(n)
*
      prdot2_8    = dble(F_xcto(n)) * dble(F_xctm_m(n)) +
     %              dble(F_ycto(n)) * dble(F_yctm_m(n)) +
     %              dble(F_zcto(n)) * dble(F_zctm_m(n))
*
      F_zcto(n) = ZERO_8
      F_ycto(n) = ZERO_8
      F_xcto(n) = ZERO_8
*
      F_xctm(n) = 2.0 * ( prx_8 * prdot2_8 ) + F_xctm(n)
      F_yctm(n) = 2.0 * ( pry_8 * prdot2_8 ) + F_yctm(n)
      F_zctm(n) = 2.0 * ( prz_8 * prdot2_8 ) + F_zctm(n)
*
      enddo
      enddo
      enddo
!$omp enddo
!$omp end parallel
*
      return
      end