!-------------------------------------- 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_fgmres -  call iterative elliptic solver 
*
#include "model_macros_f.h"

      subroutine sol_fgmres( F_w2_8, F_w1_8, iln, nil, njl, 1,7
     $                          Minx, Maxx, Miny, Maxy, Nk )
      implicit none
*
      integer iln, nil, njl, Minx, Maxx, Miny, Maxy, Nk
      real*8 F_w1_8 (Minx:Maxx,Miny:Maxy,Nk),
     $       F_w2_8 (Minx:Maxx,Miny:Maxy,Nk)
*
*author
*       Abdessamad Qaddouri - December  2006
*
*revision
* v3_30 - Qaddouri A.       - initial version
* v3_31 - Tanguay M.        - Introduce time extrapolation
*
#include "glb_pil.cdk"
#include "glb_ld.cdk"
#include "lun.cdk"
#include "sol.cdk"
#include "schm.cdk"
#include "orh.cdk"
*
      integer niloc,njloc,nloc
      integer halox, haloy, minx1, maxx1, minx2, maxx2
      real*8, dimension(:  ), allocatable :: wk11,wk22,rhs1,sol1
      real*8, dimension(:,:), allocatable :: vv_8,ww_8
      integer icode,iout,ischmi
*
*     ---------------------------------------------------------------
*
      niloc = (Nil-pil_e)-(1+pil_w)+1
      njloc = (Njl-pil_n)-(1+pil_s)+1
      nloc  = niloc*njloc*Schm_nith
*
      allocate (vv_8(nloc,sol_im+1),ww_8(nloc,sol_im),wk11(nloc),
     $          wk22(nloc),rhs1(nloc),sol1(nloc))
*
      halox = 1
      haloy = halox
*
      minx1 = 1-halox
      maxx1 = nil+halox
      minx2 = 1-haloy
      maxx2 = njl+haloy
*
      sol_i0 = 1   + pil_w
      sol_in = nil - pil_e
      sol_j0 = 1   + pil_s
      sol_jn = njl - pil_n
*
      call tab_vec ( F_w1_8 , Minx,Maxx,Miny,Maxy, Schm_nith,
     $               rhs1   , sol_i0,sol_in,sol_j0,sol_jn, +1 )
*
      sol1  = 0.
      iout  = 0
      icode = 0
*
      if ((iln.eq.Schm_itnlh).and.(Orh_icn.eq.Schm_itcn.or..not.Orh_crank_L)) iout=lun_out
*
*-----------------------------------------------------------------------
*     F G M R E S   L O O P
*-----------------------------------------------------------------------
   1    continue
*-----------------------------------------------------------------------
      call fgmres (nloc,sol_im,rhs1,sol1,ischmi,vv_8,ww_8,wk11,wk22,
     *                                sol_eps,sol_maxits,iout,icode)
      if (icode .eq. 1) then
*
         if (sol_precond_S.eq.'DIAGO')         then
            call sol_pre_diago ( wk22, wk11, Minx,Maxx,Miny,Maxy,
     $                           nil,njl, minx1,maxx1,minx2,maxx2 )
*
         else if (sol_precond_S.eq.'JACOBI')   then
            call sol_pre_jacobi   ( wk22, wk11, niloc, njloc )
*
         else if (sol_precond_S.eq.'MULTICOL') then
            call sol_pre_multicol ( wk22, wk11, niloc, njloc,
     $                              Minx, Maxx, Miny, Maxy, nil,njl,
     $                              minx1, maxx1, minx2, maxx2, Nk )
*
         else
*
            call dcopy (nloc, wk11, 1, wk22, 1)
*
         endif
*
         goto 1
*
      else
*
         if (icode.ge.2) then
*
            call  sol_matvecs ( wk22, wk11, Minx, Maxx, Miny, Maxy,
     $                            nil,njl, minx1,maxx1,minx2,maxx2 )
            goto 1
         endif
*
      endif
*
      call tab_vec ( F_w2_8 , Minx,Maxx,Miny,Maxy, Schm_nith,
     $               sol1   , sol_i0,sol_in,sol_j0,sol_jn, -1 )
*
      deallocate (wk11,wk22,rhs1,sol1,vv_8,ww_8)
*
*     ---------------------------------------------------------------
*
      return
      end