!-------------------------------------- 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 set_opr - initialize the commons containing model operators
*
#include "model_macros_f.h"
*

      subroutine set_opr 1,11
*
      implicit none
*
*author
*     michel roch - rpn - june 1993
*
*revision
* v2_00 - Desgagne/Lee      - initial MPI version (from setopr v1_03)
* v2_11 - Desgagne M.       - vertical sponge layer
* v3_00 - Lee/Qaddouri      - Lam configuration
* v3_30 - Dugas B.          - Corriger l'allocation avec hpalloc
* v3_30 - Qaddouri A.       - add eigenmodes with definite parity
*       - Qaddouri A.       - add call for setting of preconditioning variables 
*
*object
*	This subroutine initializes the commons containing the
*	matrices used by the operators of the model
*	
*arguments
*	none
*
*implicits
#include "glb_ld.cdk"
#include "glb_pil.cdk"
#include "dcst.cdk"
#include "fft.cdk"
#include "lun.cdk"
#include "geomg.cdk"
#include "schm.cdk"
#include "opr.cdk"
#include "eigv.cdk"
#include "sol.cdk"
#include "cstv.cdk"
#include "trp.cdk"
#include "prec.cdk" 
#include "ptopo.cdk"
*
*modules
**
      real*8 ZERO_8, ONE_8, HALF_8
      parameter( ZERO_8 = 0.0 )
      parameter( ONE_8  = 1.0 )
      parameter( HALF_8 = 0.5 )
*
      integer i, j, k, k0, l, dim, err,Gni,Gnj,NSTOR
      real*8  sc_8, gdx_8
      real*8, dimension(:)  ,allocatable :: wk_8, wk2_8
*
*     ---------------------------------------------------------------
*
      call set_transpose ()
*
*C    Initialize commons for fast fourier transforms
*
      call set_fft ()
*
      dim = (trp_12smax-trp_12smin+1)*(trp_22max-trp_22min+1)*G_nj
      allocate (Sol_ai_8(dim),Sol_bi_8(dim),Sol_ci_8(dim))
*
      if (Lun_out.gt.0) write(Lun_out,1000)
*
*     Initialize projection operators to ZERO_8
*
      dim = 3*G_ni
      call hpalloc(Opr_opsxp0_8_ ,dim, err,8)
      call hpalloc(Opr_opsxp2_8_ ,dim, err,8)
      dim = 3*G_nj
      call hpalloc(Opr_opsyp0_8_ ,dim, err,8)
      call hpalloc(Opr_opsyp2_8_ ,dim, err,8)
      dim = 3*G_nk
      call hpalloc(Opr_opszp0_8_ ,dim, err,8)
      call hpalloc(Opr_opszp2_8_ ,dim, err,8)
      do i = 1, 3*G_ni
         Opr_opsxp0_8 (i) = ZERO_8
         Opr_opsxp2_8 (i) = ZERO_8
      end do
      do j = 1, 3*G_nj
         Opr_opsyp0_8 (j) = ZERO_8
         Opr_opsyp2_8 (j) = ZERO_8
      end do    
      do k = 1, 3*G_nk
         Opr_opszp0_8 (k) = ZERO_8
         Opr_opszp2_8 (k) = ZERO_8
      end do
*
*C    Allocate memory for eigenvectors
*
      if ( .not. Fft_fast_L ) then
         call hpalloc( Opr_xevec_8_, G_ni*G_ni, err, 8 )
      endif
      call hpalloc( Opr_xeval_8_,      G_ni, err, 8 )
      call hpalloc( Opr_zevec_8_, G_nk*G_nk, err, 8 )
      call hpalloc( Opr_zeval_8_,      G_nk, err, 8 )
*
*     Calculate dimension for dimension without pilot region
      Gni = G_ni-Lam_pil_w-Lam_pil_e
      Gnj = G_nj-Lam_pil_s-Lam_pil_n
*
*C    Prepare projection operators
*
*        1.1 Compute the East-West operators
*
      dim = max(Gni,Gnj)

      allocate ( wk_8(dim) )

      do i = 1+Lam_pil_w, G_ni-Lam_pil_e
         Opr_opsxp0_8(G_ni+i) = (G_xg_8(i+1)-G_xg_8(i-1))*HALF_8
      end do

      do i = 1+Lam_pil_w, G_ni-Lam_pil_e
         wk_8(i-Lam_pil_w) = G_xg_8(i+1) - G_xg_8(i)
      end do
      allocate ( wk2_8(Gni*3) )
      call set_ops8 (wk2_8,wk_8,ONE_8,G_periodx,Gni,Gni,1)
      do i=1,Gni
         Opr_opsxp2_8(i+Lam_pil_w)=wk2_8(i)
         Opr_opsxp2_8(G_ni+i+Lam_pil_w)=wk2_8(Gni+i)
         Opr_opsxp2_8(G_ni*2+i+Lam_pil_w)=wk2_8(Gni*2+i)
      enddo
      deallocate ( wk2_8 )
*
*
*        1.2 Compute the North-South operators
*
      do j = 1+Lam_pil_s, G_nj-Lam_pil_n
         Opr_opsyp0_8(G_nj+j) =
     $        sin((G_yg_8(j+1)+G_yg_8(j  ))* HALF_8)- 
     $        sin((G_yg_8(j  )+G_yg_8(j-1))* HALF_8)
      end do
*
      do j = 1+Lam_pil_s, G_nj-1-Lam_pil_n
         wk_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 set_ops8(wk2_8,wk_8,ZERO_8,G_periody,Gnj,Gnj,1)
      do j=1,Gnj
         Opr_opsyp2_8(j+Lam_pil_s)=wk2_8(j)
         Opr_opsyp2_8(G_nj+j+Lam_pil_s)=wk2_8(Gnj+j)
         Opr_opsyp2_8(G_nj*2+j+Lam_pil_s)=wk2_8(Gnj*2+j)
      enddo
      deallocate ( wk2_8 )
*
      deallocate ( wk_8 )
      call set_ops8 (Opr_opszp0_8,geomg_hz_8,ONE_8,.false.,G_nk,G_nk,2)
*
*C    Compute eigenvalues and eigenvector for the generalized
*     eigenvalue problem in East-West direction
*
*
      if ( .not. Fft_fast_L ) then
*        Compute eigenvalues when no FFT are used:

         if( .not. Eigv_parity_L) then
*
            call set_poic  ( Opr_xeval_8, Opr_xevec_8 , Opr_opsxp0_8,
     $                                       Opr_opsxp2_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(Opr_evvec_8_ ,dim, err,8)
            call hpalloc(Opr_odvec_8_ ,dim, err,8)
            call  set_poic_par ( Opr_xeval_8, Opr_evvec_8, Opr_odvec_8,
     $                                         G_xg_8(1), G_ni, NSTOR )
         endif
*
      else

*
*        Compute eigenvalues when FFT are used:
*
         sc_8 = Dcst_pi_8 / dble( Gni )
         if (Lun_debug_L) print *,'FFT sc_8=',sc_8,' Gni=',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
                Opr_xeval_8(i)    = ZERO_8
             enddo
             do i=G_ni-Lam_pil_e+1,G_ni
                Opr_xeval_8(i)    = ZERO_8
             enddo
             do i = 2+Lam_pil_w, G_ni-Lam_pil_e
                Opr_xeval_8(i) = - (2*sin(float(i-Lam_pil_w-1)*sc_8/2)/gdx_8)**2
             enddo
         else
             Opr_xeval_8(1)    = ZERO_8
             Opr_xeval_8(G_ni) = - ONE_8 / ( sc_8 ** 2. )
             do i = 1, (G_ni-1)/2
                Opr_xeval_8(2*i+1) = - (sin( dble(i) * sc_8 ) / sc_8 ) **2.
                Opr_xeval_8(2*i)   =  Opr_xeval_8(2*i+1)
             end do
         endif
*
      endif
*
      call set_oprz ()
*
      if ( sol_type_S.eq.'ITERATIF' ) then
         if (Lun_out.gt.0) write (Lun_out,1001) trim(sol_precond_S)
         if ( (sol_precond_S.eq.'JACOBI') .or. 
     $        (sol_precond_S.eq.'MULTICOL') ) call set_prec
      endif
*   
*C    Initialize common block for diffusion
*
      call hzd_set ()
*   
*C    Initialize common block for vertical sponge
*
      call vspng_set ()
*
 1000 format(/,'INITIALIZATING MODEL OPERATORS    (S/R SET_OPR)',
     %       /,'=============================================')
 1001 format(/,'WILL USE FGMRES ITERATIVE SOLVER WITH ',a,
     $         ' PRECONDITIONNER'
     %       /,'=============================================')
*
*     ---------------------------------------------------------------
*
      return
      end