!-------------------------------------- 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 sol_hcr - iterative solution of horizontal Helmholtz problem
*               cgrid model with rotated poles
*
#include "model_macros_f.h"
*

      subroutine sol_hcr (F_sol_8,F_rhs_8,F_w1_8,F_w2_8,F_dg1_8,F_dg2_8, 2,5
     $                    F_dwfft_8,iln,Minx,Maxx,Miny,Maxy,Ni,Nj,Nk)
      implicit none
*
      integer iln,Minx,Maxx,Miny,Maxy,Ni,Nj,Nk
      real*8 F_sol_8 (Minx:Maxx,Miny:Maxy,Nk),
     $       F_rhs_8 (Minx:Maxx,Miny:Maxy,Nk),
     $        F_w1_8 (Minx:Maxx,Miny:Maxy,Nk),
     $        F_w2_8 (Minx:Maxx,Miny:Maxy,Nk),
     $        F_dg1_8(*),F_dg2_8(*),F_dwfft_8(*)
*
*author
*     Desgagne/Lee ( after version v1_03 of solhcr8 )
*
*revision
* v2_00 - Desgagne/Lee       - initial MPI version
* v3_00 - Desgagne & Lee     - Lam configuration
* v3_01 - Qaddouri A.        - add call sol_parite_2
* v3_10 - Corbeil & Desgagne & Lee - AIXport+Opti+OpenMP
* v3_21 - Desgagne M.        - Optimization
* v3_30 - Qaddouri A.        - add Call to iterative solver 
*object
*
*arguments
*  Name        I/O                 Description
*----------------------------------------------------------------
* F_sol_8
*----------------------------------------------------------------
*
#include "glb_ld.cdk"
#include "glb_pil.cdk"
#include "ldnh.cdk"
#include "sol.cdk"
#include "opr.cdk"
#include "eigv.cdk"
#include "ptopo.cdk"
#include "fft.cdk"
#include "cstv.cdk"
#include "schm.cdk"
#include "trp.cdk"
#include "lctl.cdk"
#include "lun.cdk"
*
      integer i, j, k, k0, offi, offj, Gni, Gnj,NSTOR,nev
      real*8 tsolver, tsolverflop, tsolver1, tsolverflop1,
     $       abpt((maxy-miny+1)*(trp_12smax-trp_12smin+1)*G_ni)
      real*8, dimension(:),allocatable :: wk_evec_8
      real*8 con(G_nk)
*
*
*     ---------------------------------------------------------------
*
* Compute length of working vector without pilot region
*
      Gni=G_ni-Lam_pil_w-Lam_pil_e
      Gnj=G_nj-Lam_pil_s-Lam_pil_n
*
      if (.not.Fft_fast_L) then
          allocate ( wk_evec_8(Gni*Gni) )
          do j=1,Gni
          do i=1,Gni
            wk_evec_8((j-1)*Gni+i)=
     $              Opr_xevec_8((j+Lam_pil_w-1)*G_ni+i+Lam_pil_w)
          enddo
          enddo
      endif
*
!$omp parallel shared ( offi,offj,g_nk )
      offi = Ptopo_gindx(1,Ptopo_myproc+1)-1
      offj = Ptopo_gindx(3,Ptopo_myproc+1)-1
*
!$omp do
      do j=1+pil_s,Nj-pil_n
         call dgemm('N','N', (ni-pil_w-pil_e), G_nk, G_nk, 1.0D0, 
     $              F_rhs_8(1+pil_w,j,1), (Maxy-Miny+1)*(Maxx-Minx+1),
     $              Opr_zevec_8,g_nk,0.0d0,
     $              F_w1_8 (1+pil_w,j,1), (Maxy-Miny+1)*(Maxx-Minx+1))
*
         do k=1,Schm_nith
            do i = 1+pil_w, Ni-pil_e
               F_w1_8(i,j,k)= Opr_opsxp0_8(G_ni+offi+i) *
     $                        Opr_opsyp0_8(G_nj+offj+j) * F_w1_8(i,j,k) 
            enddo
         end do
      end do
!$omp enddo
*
!$omp end parallel 
*     
      if ( sol_type_S.eq.'ITERATIF' ) then

         call  sol_fgmres ( F_w2_8, F_w1_8, iln, l_ni, l_nj,
     $                           Minx, Maxx, Miny, Maxy, Nk )

      else

         if (Fft_fast_L) then
         if (G_lam) then
            
             call sol_fft8_lam ( F_w2_8, F_w1_8, Fft_pri_8,
     $                     Minx, Maxx, Miny, Maxy, ldnh_nj,
     $            trp_12smin, trp_12smax, Schm_nith, trp_12sn , 
     $            G_ni, G_nj, trp_22min , trp_22max, trp_22n  ,
     $            trp_12smin, trp_12smax, trp_22min, trp_22max,G_nj,
     $            Ptopo_npex, Ptopo_npey, Sol_ai_8,Sol_bi_8,Sol_ci_8,
     $            F_dg2_8,F_dwfft_8)
         else
             call sol_fft8_2 ( F_w2_8, F_w1_8, Fft_pri_8,
     $                     Minx, Maxx, Miny, Maxy, ldnh_nj,
     $            trp_12smin, trp_12smax, Schm_nith, trp_12sn , 
     $            G_ni, G_nj, trp_22min , trp_22max, trp_22n  ,
     $            trp_12smin, trp_12smax, trp_22min, trp_22max,G_nj,
     $            Ptopo_npex, Ptopo_npey, Sol_ai_8,Sol_bi_8,Sol_ci_8,
     $            F_dg2_8,F_dwfft_8)
         endif
*
         else
*
         if(.not. Eigv_parity_L) then
*
            call sol_mxma8_2 ( F_w2_8, F_w1_8, wk_evec_8, 
     $                      Minx, Maxx, Miny, Maxy,ldnh_nj,
     $            trp_12smin, trp_12smax, Schm_nith, trp_12sn , 
     $            G_ni, G_nj, trp_22min , trp_22max, trp_22n  ,
     $            trp_12smin, trp_12smax, trp_22min, trp_22max,G_nj,
     $            Ptopo_npex, Ptopo_npey, Sol_ai_8,Sol_bi_8,Sol_ci_8,
     $            F_dg1_8,F_dg2_8,F_dwfft_8)
         else
*
            nev= (G_ni+2)/2
            NSTOR = nev + ( 1 - mod(nev,2) )
*
            call sol_parite_2( F_w2_8, F_w1_8, Opr_evvec_8, Opr_odvec_8,
     $                       Minx, Maxx, Miny, Maxy, l_nj,
     $            trp_12smin, trp_12smax, Schm_nith, trp_12sn ,
     $            G_ni, G_nj, trp_22min , trp_22max, trp_22n  ,
     $            trp_12smin, trp_12smax, trp_22min, trp_22max,G_nj,
     $            Ptopo_npex, Ptopo_npey, Sol_ai_8,Sol_bi_8,Sol_ci_8,
     $            F_dg1_8,F_dg2_8,F_dwfft_8,Abpt,NSTOR,nev)
*
         endif
         endif
*
      endif
*
*     inverse projection 
*
      do k=1,G_nk-1
         con(k) = 1.
      enddo
      con(G_nk) = -1./Cstv_hco0_8
*
!$omp parallel shared ( g_nk )
!$omp do
      do j=1+pil_s,Nj-pil_n
         F_w2_8(1+pil_w:Ni-pil_e,j,G_nk) =  
     $                  con(G_nk)*F_w1_8(1+pil_w:Ni-pil_e,j,G_nk)
         call dgemm('N','T', (ni-pil_w-pil_e), G_nk, G_nk, 1.0D0, 
     $              F_w2_8 (1+pil_w,j,1), (Maxy-Miny+1)*(Maxx-Minx+1),
     $              Opr_zevec_8,g_nk,0.0d0,
     $              F_sol_8(1+pil_w,j,1), (Maxy-Miny+1)*(Maxx-Minx+1))
      enddo
!$omp enddo
!$omp end parallel
*
      if (.not. Fft_fast_L) deallocate (wk_evec_8)
*
*     ---------------------------------------------------------------
* 
      return
      end