!-------------------------------------- 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 hzd_set - Compute diffusion operator matrix on G, U and V grids
*
#include "model_macros_f.h"
*

      subroutine hzd_set 1,14
      implicit none
*
*author    
*     J.P. Toviessi - CMC - Jan 1999
*
*revision
* v2_00 - Desgagne M.       - initial MPI version
* v2_10 - Qaddouri&Desgagne - higher order diffusion operator
* v2_11 - Desgagne M.       - remove vertical modulation
* v2_31 - Qaddouri A.       - remove stkmemw and correction to yp2
* v3_00 - Qaddouri & Lee    - Lam configuration
* v3_01 - Toviessi J. P.    - add eigenmodes with definite parity
* v3_01 - Lee V.            - add setup for horizontal sponge
* v3_20 - Qaddouri/Toviessi - variable higher order diffusion operator
* v3_20 - Tanguay M.        - 1d higher order diffusion operator
* v3_30 - Tanguay M.        - activate Hzd_type_S='EXPLICIT' 
*
*object
*     see id section
*       
*arguments
*     none
*
*implicits
#include "glb_ld.cdk"
#include "glb_pil.cdk"
#include "geomg.cdk"
#include "hzd.cdk"
#include "eigv.cdk"
#include "dcst.cdk"
#include "fft.cdk"
#include "trp.cdk"
#include "opr.cdk"
#include "cstv.cdk"
#include "lun.cdk"
#include "vspng.cdk"
#include "hspng.cdk"
#include "ptopo.cdk"
*      
      character*11 mess
      integer i,j,k,istat,dim,dpwr,Gni,Gnj,NSTOR,dimx,dimy
      real*8, dimension(:) , allocatable :: wk1_8,wk2_8,
     $                                      h0_8,h0t_8,h2_8,h2t_8
      real*8 c_8, gdx_8
      real*8 ZERO_8, ONE_8, HALF_8
      PARAMETER( ZERO_8 = 0.0 , ONE_8  = 1.0 , HALF_8 = 0.5 )
**
*     ---------------------------------------------------------------
*
*     Set the length of the work matrices to the region without 
*     the pilot region for LAM
*
      Gni = G_ni-Lam_pil_w-Lam_pil_e
      Gnj = G_nj-Lam_pil_s-Lam_pil_n
*
      Hzd_fact_L = .false.
      Hzd_ho_L   = .false.
      mess = 'NIL'
*
      if ((Hzd_type_S.eq."FACT").and.
     $   ((Cstv_uvdf_8.gt.0.).or.(Cstv_phidf_8.gt.0)) ) then
         Hzd_fact_L = .true.
         mess       = 'FACTORIZED'
      endif
*
      if ((Hzd_lnr.gt.0.).and. 
     $    ((Hzd_type_S.eq."HO").or.(Hzd_type_S.eq."EXPLICIT"))) then
         Hzd_ho_L   = .true.
         mess       = 'HIGH ORDER'
      endif
*
      if (Lun_out.gt.0) write(Lun_out,1001) mess
*
      if ((Hzd_fact_L.or.Hzd_ho_L).or.(Vspng_nk.gt.0)
     $                            .or.(Hspng_nj.gt.0) ) then
*
         dimx = 3*G_ni*2
         dimy = 3*G_nj*2
         call hpalloc(VAR_PTR(Hzd_xp0_8)      ,dimx,istat,1)
         call hpalloc(VAR_PTR(Hzd_xp2_8)      ,dimx,istat,1)
         call hpalloc(VAR_PTR(Hzd_yp0_8)      ,dimy,istat,1)
         call hpalloc(VAR_PTR(Hzd_yp2_8)      ,dimy,istat,1)
         call hpalloc(VAR_PTR(Hzd_yp2su_8 )   ,dimy,istat,1)
*
         if (Hzd_1d_L) then
            call hpalloc(VAR_PTR(Hz1d_yp2_8   ),dimy,istat,1) 
            call hpalloc(VAR_PTR(Hz1d_yp2su_8 ),dimy,istat,1)
         endif
*
         allocate (h0_8(G_nj*2),h0t_8(G_nj*2),h2_8(G_nj*2)
     %                                      ,h2t_8(G_nj*2))
*
*______________________________________________________________________
* NOTES relative a la difusion hizon. variable:
*------
*  1) calcul de xp0 et xp2
*     Comme g0_8 == g0t_8 == 1.0 on a alors Hzd_xp0sv_8 == opr_opsxp0
*                               Hzd_xp2sv_8 == opr_opsxp2
*        et                     hzd_xp0_8 ne change pas  pour U
*                               Hzd_xp2_8  == hzd_xp2_8 pour U
*        donc on ne reacule pas hzd_xp0_8
*        et on utilise les valeurs calculees de opr_opsxp0, 
*        opr_opsxp2,hzd_xp2_8
*        on utilise aussi les valeurs Hzd_xeval_8 qui auraient ete
*        Hzd_xevalsv_8 et Hzd_xeval_8 si g0_8 # 1.0 # g0t_8.
*        de plus on utilise aussi les valeurs Hzd_xevec_8 qui auraient 
*        ete Hzd_xevecsv_8 et Hzd_xevec_8 si g0_8 # 1.0 # g0t_8.
*        Hzd_xevecsv_8,Hzd_xevec_8
*  2) calcul de yp2 (remarque yp0 ne change pas ): On doit recalculer
*      opr_opsyp2 dans Hzd_yp2su ( pour S et U) et Hzd_yp2 pour V
*
         if (Hzd_difva_L) then
           do i = 1, G_nj
             h0_8(i)  = ONE_8 / (ONE_8 + 2*((cos(G_yg_8(i)))**2))
             h0t_8(i) = ONE_8 / (ONE_8 + 2*((cos(G_yg_8(i)))**2))
*            h0 == h2 and h2t == h0t
             h2_8(i)  = h0_8(i)
             h2t_8(i) = h0t_8(i)
           enddo

           do j = 1, 3*G_nj
              Hzd_yp2su_8 (j) = ZERO_8
           end do
*
         else
           do i = 1, G_nj
             h0_8(i)  = ONE_8
             h0t_8(i) = ONE_8
*            h0 == h2 and h2t == h0t
             h2_8(i)  = h0_8(i)
             h2t_8(i) = h0t_8(i)
           enddo
         endif
*_______________________________________________________________________

         do i = 1, 3*G_ni
            Hzd_xp0_8(i) = ZERO_8
         end do
         do i = 1+Lam_pil_w, G_ni-Lam_pil_e
            Hzd_xp0_8(G_ni+i) = G_xg_8(i+1) - G_xg_8(i)
         end do

         allocate ( wk1_8(Gni)   )
         allocate ( wk2_8(Gni*3) )

         do i = 1, 3*G_ni
            Hzd_xp2_8(i) = ZERO_8
         end do

         do i = 1+Lam_pil_w, G_ni-Lam_pil_e
            wk1_8(i-Lam_pil_w) = (G_xg_8(i+2)-G_xg_8(i)) * HALF_8
         end do
         call set_ops8 (wk2_8,wk1_8,ONE_8,G_periodx,Gni,Gni,1)
*
         do i=1,Gni
            Hzd_xp2_8(i+Lam_pil_w)=wk2_8(i)
            Hzd_xp2_8(G_ni+i+Lam_pil_w)=wk2_8(Gni+i)
            Hzd_xp2_8(G_ni*2+i+Lam_pil_w)=wk2_8(Gni*2+i)
         enddo
*
         deallocate (wk1_8,wk2_8)
*
         do j = 1, 3*G_nj
            Hzd_yp0_8(j) = ZERO_8
            Hzd_yp2_8(j) = ZERO_8
         end do    
         do j = 1+Lam_pil_s, G_nj-Lam_pil_n
            Hzd_yp0_8(G_nj+j) = sin(G_yg_8(j+1))-sin(G_yg_8(j))
         end do
*
         j=1+Lam_pil_s
         Hzd_yp2_8(2*G_nj+j)= ((cos(G_yg_8(j+1))**2)*h2_8(j+1))/(
     %                      sin((G_yg_8(j+2)+G_yg_8(j+1  ))* HALF_8)-
     %                      sin((G_yg_8(j  )+G_yg_8(j+1))* HALF_8))
         Hzd_yp2_8(G_nj+j) =-Hzd_yp2_8(2*G_nj+j)

         do j = 2+Lam_pil_s, G_njv-1-Lam_pil_n
            Hzd_yp2_8(2*G_nj+j)= ((cos(G_yg_8(j+1))**2)*h2_8(j+1))/(
     %                          sin((G_yg_8(j+2)+G_yg_8(j+1))* HALF_8)-
     %                          sin((G_yg_8(j+1)+G_yg_8(j  ))* HALF_8))
            Hzd_yp2_8(j) = ((cos(G_yg_8(j))**2)*h2_8(j-1))/(
     %                      sin((G_yg_8(j+1)+G_yg_8(j  ))* HALF_8)-
     %                      sin((G_yg_8(j  )+G_yg_8(j-1))* HALF_8))
            Hzd_yp2_8(G_nj+j) = - (Hzd_yp2_8(j) + Hzd_yp2_8(2*G_nj+j))
         enddo

         j=G_njv-Lam_pil_n
         Hzd_yp2_8(j) = h2_8(j-1)*Hzd_yp2_8(2*G_nj+j-1)
         Hzd_yp2_8(G_nj+j) = - (Hzd_yp2_8(j) + Hzd_yp2_8(2*G_nj+j))
*
         if (Hzd_1d_L) then
            do j=1,3*G_nj
               Hz1d_yp2_8(j) = ZERO_8
            enddo
         endif
*
*______________________________________________________________________
*
         if (Hzd_difva_L) then
            allocate ( wk1_8(Gnj) )
            do j = 1+Lam_pil_s, G_nj-1-Lam_pil_n
               wk1_8(j-Lam_pil_s) = (sin (G_yg_8(j+1))-sin(G_yg_8(j))) /
     $             (cos ((G_yg_8(j+1)+G_yg_8(j))*HALF_8)**2)
            end do
            allocate ( wk2_8(Gnj*3) )
            call hzd_set_ops8(wk2_8,wk1_8,ZERO_8,G_periody,Gnj,Gnj,1,h2t_8)
            do j=1,Gnj
               Hzd_yp2su_8(j+Lam_pil_s)=wk2_8(j)
               Hzd_yp2su_8(G_nj+j+Lam_pil_s)=wk2_8(Gnj+j)
               Hzd_yp2su_8(G_nj*2+j+Lam_pil_s)=wk2_8(Gnj*2+j)
            enddo
            deallocate (wk1_8, wk2_8 )
*
            if (Hzd_1d_L) then
                do j=1,3*G_nj
                   Hz1d_yp2su_8(j) = ZERO_8
                enddo
            endif
*
         endif
*______________________________________________________________________
*
         if (Hzd_fact_L.or.Hspng_nj.gt.0) then
*
            if (Lun_out.gt.0.and.Hspng_nj.gt.0) write(Lun_out,1004)
*
            call hpalloc(VAR_PTR(Hzd_opsxp0_8)   ,  G_ni*2,istat,1)
            call hpalloc(VAR_PTR(Hzd_opsyp0_8)   ,  G_nj*2,istat,1)
*
            do i = 1,G_ni
               Hzd_opsxp0_8 (i) = ZERO_8
            enddo
            do i = 1+Lam_pil_w, G_ni-Lam_pil_e
               Hzd_opsxp0_8 (i) = G_xg_8(i+1) - G_xg_8(i)
            end do
*
            do j = 1,G_nj
               Hzd_opsyp0_8 (j) = ZERO_8
            enddo
            do j = 1+Lam_pil_s, G_nj-Lam_pil_n
               Hzd_opsyp0_8 (j) = sin(G_yg_8(j+1))-sin(G_yg_8(j))
            end do
*
         endif
*
*     Compute eigenvalues and eigenvector for high-order diffusion.
*              Eigenvalue problem in East-West direction
*        -------------------------------------------------------
         if (Hzd_ho_L) then
*
            call hpalloc( Hzd_xeval_8_, G_ni*2, istat, 1 )
            if ( .not. Fft_fast_L ) then

              if( .not. Eigv_parity_L) then
                 call hpalloc( Hzd_xevec_8_, G_ni*G_ni*2, istat, 1 )
                 call set_poic  ( Hzd_xeval_8, Hzd_xevec_8 , Hzd_xp0_8,
     $                                            Hzd_xp2_8, Gni, G_ni )
              else
*                Eigenmodes with definite parity
                 NSTOR = (G_ni+2)/2 + ( 1 - mod((G_ni+2)/2,2) )
                 dim = NSTOR*NSTOR
                 call hpalloc( Hzd_evvec_8_, dim, istat, 8 )
                 call hpalloc( Hzd_odvec_8_, dim, istat, 8 )
                 call set_poic_par_U (Hzd_xeval_8, Hzd_evvec_8 ,
     $                              Hzd_odvec_8, G_xg_8(1), G_ni, NSTOR)
              endif
            else
               c_8 = Dcst_pi_8 / dble( Gni )
               if (G_lam) then
                   gdx_8 = (G_xg_8(G_ni-Lam_pil_e)-G_xg_8(Lam_pil_w) )/dble(Gni)
                   if(Lun_debug_L)print *,'gdx=',gdx_8
                   do i=1,1+Lam_pil_w
                      Hzd_xeval_8(i)    = ZERO_8
                   enddo
                   do i=G_ni-Lam_pil_e+1,G_ni
                      Hzd_xeval_8(i)    = ZERO_8
                   enddo
                   do i = 2+Lam_pil_w, G_ni-Lam_pil_e
                      Hzd_xeval_8(i) = 
     $                    - (2*sin(float(i-Lam_pil_w-1)*c_8/2)/gdx_8)**2
                   enddo
               else
                   Hzd_xeval_8(1)    =   ZERO_8
                   Hzd_xeval_8(G_ni) = - ONE_8 / ( c_8 ** 2. )
                   do i = 1, (G_ni-1)/2
                      Hzd_xeval_8(2*i+1) = - (sin(dble(i) * c_8) / c_8) **2.
                      Hzd_xeval_8(2*i)   =  Hzd_xeval_8(2*i+1)
                   end do
               endif
            endif
*
*     initialize operator nord_south for U, V and scalar grids
*                              for high-order diffusion-solver
*
            dpwr= Hzd_pwr / 2
            dim = (trp_22max-trp_22min+1)
            dim = dim * G_nj * dpwr*dpwr
            call hpalloc(Hzd_au_8_    , dim, istat, 8)
            call hpalloc(Hzd_cu_8_    , dim, istat, 8)
            call hpalloc(Hzd_deltau_8_, dim, istat, 8)
            call hpalloc(Hzd_av_8_    , dim, istat, 8)
            call hpalloc(Hzd_cv_8_    , dim, istat, 8)
            call hpalloc(Hzd_deltav_8_, dim, istat, 8)
            call hpalloc(Hzd_as_8_    , dim, istat, 8)
            call hpalloc(Hzd_cs_8_    , dim, istat, 8)
            call hpalloc(Hzd_deltas_8_, dim, istat, 8)
*
            if (Hzd_1d_L) then
                dimx= G_ni * G_nj * dpwr*dpwr
                call hpalloc(Hz1d_deltau_8_, dimx, istat, 8)
                call hpalloc(Hz1d_deltav_8_, dimx, istat, 8)
                call hpalloc(Hz1d_deltas_8_, dimx, istat, 8)
            endif
*
            c_8= 1.E+32
            if(max(G_ni,2*G_nj).eq.G_ni) then
               do i = 1+Lam_pil_w, G_ni-Lam_pil_e
                  c_8 = min(c_8,(G_xg_8(i+1) - G_xg_8(i)))
               end do
            else
               do j= 1+Lam_pil_s,G_nj-Lam_pil_n
                  c_8 = min(c_8,(G_yg_8(j+1) - G_yg_8(j)))
               end do
            endif
            hzd_lnR   = log(1.- Hzd_lnR)
            Hzd_cdiff = (2./c_8)**Hzd_pwr / (-Hzd_lnR)
            if (Lun_out.gt.0) then
               write(Lun_out,1010) 
     $         (Dcst_rayt_8**Hzd_pwr)/(Cstv_dt_8*Hzd_cdiff),Hzd_pwr
            endif
            allocate (wk1_8(3*G_nj),wk2_8(3*G_nj))
            do j = 1,3*G_nj
               wk1_8 (j) = ZERO_8
               wk2_8 (j) = ZERO_8
            enddo
            do j = 1+Lam_pil_s, G_nj-Lam_pil_n
               wk1_8 (       j) = ZERO_8
               wk1_8 (  G_nj+j) = Opr_opsyp0_8(G_nj+j)*h0_8(j)/cos(G_yg_8(j))**2
               wk1_8 (2*G_nj+j) = ZERO_8
               wk2_8 (       j) = ZERO_8
               wk2_8 (  G_nj+j) = (sin(G_yg_8(j+1))-sin(G_yg_8(j)))*h0t_8(j)/
     %               ( cos ( (G_yg_8(j+1  )+G_yg_8(j))* HALF_8) **2)
               wk2_8 (2*G_nj+j) = ZERO_8
            end do
*
            if (Lun_debug_L) print *,'HZD_DELPWR on V'
            call hzd_delpwr (Hzd_av_8,Hzd_cv_8,Hzd_deltav_8,dpwr,
     $                       trp_22min,trp_22max,G_nj,trp_22n,trp_22n0,
     $                       Hzd_yp0_8,Hzd_yp2_8,wk2_8,
     $                       Opr_xeval_8, Hzd_cdiff)
*
            if (Hzd_difva_L) then
               if (Lun_debug_L) print *,'HZD_DELPWR on U'
               call hzd_delpwr (Hzd_au_8,Hzd_cu_8,Hzd_deltau_8,dpwr,
     $                       trp_22min,trp_22max,G_nj,trp_22n,trp_22n0,
     $                       Opr_opsyp0_8,Hzd_yp2su_8,wk1_8,
     $                       Hzd_xeval_8, Hzd_cdiff)
*
               if (Lun_debug_L) print *,'HZD_DELPWR on scalar'
               call hzd_delpwr (Hzd_as_8,Hzd_cs_8,Hzd_deltas_8,dpwr,
     $                       trp_22min,trp_22max,G_nj,trp_22n,trp_22n0,
     $                       Opr_opsyp0_8,Hzd_yp2su_8,wk1_8,
     $                       Opr_xeval_8, Hzd_cdiff)
            else
               if (Lun_debug_L) print *,'HZD_DELPWR on U'
               call hzd_delpwr (Hzd_au_8,Hzd_cu_8,Hzd_deltau_8,dpwr,
     $                       trp_22min,trp_22max,G_nj,trp_22n,trp_22n0,
     $                       Opr_opsyp0_8,Opr_opsyp2_8,wk1_8,
     $                       Hzd_xeval_8, Hzd_cdiff)
*
               if (Lun_debug_L) print *,'HZD_DELPWR on scalar'
               call hzd_delpwr (Hzd_as_8,Hzd_cs_8,Hzd_deltas_8,dpwr,
     $                       trp_22min,trp_22max,G_nj,trp_22n,trp_22n0,
     $                       Opr_opsyp0_8,Opr_opsyp2_8,wk1_8,
     $                       Opr_xeval_8, Hzd_cdiff)
            endif
*
            if (Hzd_1d_L) then
*
               if (Lun_debug_L) print *,'HZD_DELPWR 1D on V'
               call hzd_delpwr_1d (Hz1d_deltav_8,dpwr,G_ni,G_nj,
     $                          Hzd_yp0_8,Hz1d_yp2_8,wk2_8,
     $                          Opr_xeval_8, Hzd_cdiff)
*
               if (Hzd_difva_L) then
                  if (Lun_debug_L) print *,'HZD_DELPWR 1D on U'
                  call hzd_delpwr_1d (Hz1d_deltau_8,dpwr,G_ni,G_nj,
     $                          Opr_opsyp0_8,Hz1d_yp2su_8,wk1_8,
     $                          Hzd_xeval_8, Hzd_cdiff)
*
                  if (Lun_debug_L) print *,'HZD_DELPWR 1D on scalar'
                  call hzd_delpwr_1d (Hz1d_deltas_8,dpwr,G_ni,G_nj,
     $                          Opr_opsyp0_8,Hz1d_yp2su_8,wk1_8,
     $                          Opr_xeval_8, Hzd_cdiff)
               else
                  if (Lun_debug_L) print *,'HZD_DELPWR 1D on U'
                  call hzd_delpwr_1d (Hz1d_deltau_8,dpwr,G_ni,G_nj,
     $                          Opr_opsyp0_8,Opr_opsyp2_8,wk1_8,
     $                          Hzd_xeval_8, Hzd_cdiff)
*
                  if (Lun_debug_L) print *,'HZD_DELPWR 1D on scalar'
                  call hzd_delpwr_1d (Hz1d_deltas_8,dpwr,G_ni,G_nj,
     $                          Opr_opsyp0_8,Opr_opsyp2_8,wk1_8,
     $                          Opr_xeval_8, Hzd_cdiff)
               endif
*
            endif
*
            deallocate (wk1_8,wk2_8)
         endif
         deallocate (h0_8,h0t_8,h2_8,h2t_8)
      else
         if (Lun_out.gt.0) write(Lun_out,1003)
      endif
*
 1001 format(/,'INITIALIZATING ',a,' HORIZONTAL DIFFUSION ', 
     $         '(S/R HZD_SET)',/,60('='))
 1003 format(/,'NO HORIZONTAL DIFFUSION REQUIRED',/,32('='))
 1004 format(/,'INITIALIZATING FACTORIZED HORIZONTAL SPONGE', 
     $         '(S/R HZD_SET)',/,60('='))
 1010 format (3X,'Diffusion Coefficient = ',e22.14,' m**',i1,'/sec' )
*
*     ---------------------------------------------------------------
      return
      end