!-------------------------------------- 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 --------------------------------------
!

      subroutine sveigen(n,nev,ncv,Eval,v,tol) 1,8
#if defined (DOC)
*
***s/r sveigen  - Driver for ARPACK solver to calculate singular vectors
*
*Author  : M. Buehner July, 2002
c-----------------------------------------------------------------------
c SUBROUTINE: SVEIGEN
c Driver subroutine for calling Lanczos algorithm to calculate
c Singular Vectors (BGSV,ANLSV, or TESV)
c
c\Routines called:
c     ssaupd  ARPACK reverse communication interface routine.
c     sseupd  ARPACK routine that returns Ritz values and (optionally)
c             Ritz vectors.
c     snrm2   Level 1 BLAS that computes the norm of a vector.
c     saxpy   Level 1 BLAS that computes y <- alpha*x+y.
c     avhess  Matrix vector multiplication routine that computes A*x.
c     tv      Matrix vector multiplication routine that computes T*x,
c             where T is a tridiagonal matrix.  It is used in routine
c             avhess.
c n   - dimension of space
c nev - number of modes to compute
c ncv - twice the number of modes ??
c-----------------------------------------------------------------------
#endif
c
c Arguments
      IMPLICIT NONE
#include "comdim.cdk"
#include "comlun.cdk"
#include "comcva.cdk"
#include "comvfiles.cdk"
#include "comsv.cdk"
      Integer n,nev,ncv,ireslun,fnom,fclos
      Real*8 v(n,ncv), Eval(nev)
      Real*8 tol
c
c     | Local Arrays |
      Real*8 workl(ncv*(ncv+8)),workd(3*n),
     &       d(ncv,2), resid(n),ax(n)
      logical          select(ncv)
c
      POINTER (PTworkl,workl)
      POINTER (PTworkd,workd)
      POINTER (PTd,d)
      POINTER (PTresid,resid)
      POINTER (PTax,ax)
      POINTER (PTselect,select)
c
      integer          iparam(11), ipntr(11),ilen
c
c     | Local Scalars |
      character        bmat*1, which*2
      integer          ido, lworkl, info, ierr, j,
     &                 nconv, maxitr, mode, ishfts,remove_c
      logical          rvec
      Real*8
     &                 sigma
c     | Parameters |
      Real*8
     &                 zero
      parameter        (zero = 0.0E+0)
c     | BLAS & LAPACK routines used |
      Real*8
     &                 dnrm2
      external         dnrm2, daxpy, remove_c, rs_dssave
c     | Intrinsic function |
      intrinsic        dabs
c
      write(nulout,*) 'ENTERING SVEIGEN'
      call vflush(nulout)
c
c Allocate local arrays
c
      CALL HPCHECK(IERR)
      IF(IERR.EQ.0)THEN
        WRITE(NULOUT,FMT='(5X,10("*")," No memory overflow"
     S        ," detected ",10("*"))')
      ELSE
        write(NULOUT,*) 'PROBLEM!!!!!'
      ENDIF
      ilen = ncv*(ncv+8)
      CALL HPALLOC(PTworkl  ,MAX(ILEN,1),IERR,8)
      ilen = 3*n
      CALL HPALLOC(PTworkd  ,MAX(ILEN,1),IERR,8)
      ilen = 2*ncv
      CALL HPALLOC(PTd      ,MAX(ILEN,1),IERR,8)
      ilen = n
      CALL HPALLOC(PTresid  ,MAX(ILEN,1),IERR,8)
      ilen = n
      CALL HPALLOC(PTax     ,MAX(ILEN,1),IERR,8)
      ilen = ncv
      CALL HPALLOC(PTselect ,MAX(ILEN,1),IERR,8)
c
      CALL HPCHECK(IERR)
      IF(IERR.EQ.0)THEN
        WRITE(NULOUT,FMT='(5X,10("*")," No memory overflow"
     S        ," detected ",10("*"))')
      ELSE
        write(NULOUT,*) 'PROBLEM!!!!!'
      ENDIF
      write(nulout,*) 'STARTING SVeigen'
      bmat = 'I'
      which = 'LM'
c
      lworkl = ncv*(ncv+8)
c      tol = zero
      info = 0
      ido = 0
c
      ishfts = 1
      maxitr = 300
      mode   = 1
c
      iparam(1) = ishfts
      iparam(3) = maxitr
      iparam(7) = mode
c
c     %-------------------------------------------%
c     | M A I N   L O O P (Reverse communication) |
c     %-------------------------------------------%
c
 10   continue
c
c Manage a possible restart of the Lanczos algorithm
c
      if(niterjob.gt.0) then
c
        write(nulout,*) 'RESTART: iteration ',nsim3d,' out of ',niterjob
        write(nulout,*) 'JOBCODE= ',njobcode
        call vflush(nulout)
c
        if(mod(nsim3d,niterjob).eq.0.and..not.LRESTART.and.nsim3d.gt.0) then
c write restart file (internal contents origininally in SAVE statements plus main work arrays)
          ireslun=0
          ierr = fnom(ireslun,crestart, 'FTN+SEQ+UNF' , 0)
          njobcode=1
          WRITE(ireslun) nsim3d,njobcode,ido,resid,v,iparam,ipntr,workd,workl,info
          call rs_dssave(.false.,ireslun)
          ierr = fclos(ireslun)
          lrestart=.true.
          WRITE(NULOUT,*) 'GOING TO RESTART!!!'
          call vflush(nulout)
          return
        endif
c
        if(LRESTART.and.NJOBCODE.eq.1) then
c read restart file
          ireslun=0
          ierr = fnom(ireslun,crestart,'FTN+SEQ+UNF+OLD+R/O',0)
          READ(ireslun) nsim3d,njobcode,ido,resid,v,iparam,ipntr,workd,workl,info
          call rs_dssave(.true.,ireslun)
          ierr = fclos(ireslun)
          ierr = remove_c(crestart)
          write(nulout,*) 'READ RESTART FILE: NSIM,JOBCODE=',nsim3d,njobcode
          call vflush(nulout)
          lrestart=.false.
        endif
c
      endif
c
c        %---------------------------------------------%
c        | Repeatedly call the routine SSAUPD and take |
c        | actions indicated by parameter IDO until    |
c        | either convergence is indicated or maxitr   |
c        | has been exceeded.                          |
c        %---------------------------------------------%
c
c        write(nulout,*)'CALLING SSAUPD'
         call dsaupd ( ido, bmat, n, which, nev, tol, resid,
     &                 ncv, v, n, iparam, ipntr, workd, workl,
     &                 lworkl, info )
c
         if (ido .eq. -1 .or. ido .eq. 1) then
c
            call avsv (n, workd(ipntr(1)), workd(ipntr(2)))
c
            go to 10
c
         end if
c
c     write(nulout,*)'finished looping on ssaupd'
      if ( info .lt. 0 ) then
c
         write(nulout,*)  ' '
         write(nulout,*)  ' Error with _saupd, info = ', info
         write(nulout,*)  ' Check documentation in _saupd '
         write(nulout,*)  ' '
c
      else
c
         rvec = .true.
c
         call dseupd ( rvec, 'All', select, d, v, n, sigma,
     &        bmat, n, which, nev, tol, resid, ncv, v, n,
     &        iparam, ipntr, workd, workl, lworkl, ierr )
c
c        | Eigenvalues are returned in the first column |
c        | of the two dimensional array D and the       |
c        | corresponding eigenvectors are returned in   |
c        | the first NEV columns of the two dimensional |
c        | array V if requested.                        |
c
         if ( ierr .ne. 0) then
c
             write(nulout,*)  ' '
             write(nulout,*)  ' Error with _seupd, info = ', ierr
             write(nulout,*)  ' Check the documentation of _seupd. '
             write(nulout,*)  ' '
c
         else
c
             nconv =  iparam(5)
             do 20 j=1, nconv
c
c               %---------------------------%
c               | Compute the residual norm |
c               |                           |
c               |   ||  A*x - lambda*x ||   |
c               |                           |
c               | for the NCONV accurately  |
c               | computed eigenvalues and  |
c               | eigenvectors.  (iparam(5) |
c               | indicates how many are    |
c               | accurate to the requested |
c               | tolerance)                |
c               %---------------------------%
c
                call avsv(n, v(1,j), ax)
                call daxpy(n, -d(j,1), v(1,j), 1, ax, 1)
                d(j,2) = dnrm2(n, ax, 1)
                d(j,2) = d(j,2) / dabs(d(j,1))
c
 20          continue
c
             call dmout(6, nconv, 2, d, ncv, -6,
     &            'Ritz values and relative residuals')
         end if
c
         if ( info .eq. 1) then
            write(nulout,*)  ' '
            write(nulout,*)  ' Maximum number of iterations reached.'
            write(nulout,*)  ' '
         else if ( info .eq. 3) then
            write(nulout,*)  ' '
            write(nulout,*)  ' No shifts could be applied during implicit',
     &               ' Arnoldi update, try increasing NCV.'
            write(nulout,*)  ' '
         end if
c
         write(nulout,*)  ' '
         write(nulout,*)  ' _SDRV1 '
         write(nulout,*)  ' ====== '
         write(nulout,*)  ' '
         write(nulout,*)  ' Size of the matrix is ', n
         write(nulout,*)  ' The number of Ritz values requested is ', nev
         write(nulout,*)  ' The number of Arnoldi vectors generated',
     &            ' (NCV) is ', ncv
         write(nulout,*)  ' What portion of the spectrum: ', which
         write(nulout,*)  ' The number of converged Ritz values is ',
     &              nconv
         write(nulout,*)  ' The number of Implicit Arnoldi update',
     &            ' iterations taken is ', iparam(3)
         write(nulout,*)  ' The number of OP*x is ', iparam(9)
         write(nulout,*)  ' The convergence criterion is ', tol
         write(nulout,*)  ' '

         do j=1,nev
           Eval(j)=d(j,1)
         enddo
c
      end if
c
 9000 continue
c
c Deallocate local arrays
c
      CALL HPCHECK(IERR)
      IF(IERR.EQ.0)THEN
        WRITE(NULOUT,FMT='(5X,10("*")," No memory overflow"
     S        ," detected ",10("*"))')
      ELSE
        write(NULOUT,*) 'PROBLEM!!!!!'
      ENDIF
      write(nulout,*) 'before dealloc'
      call vflush(nulout)
      CALL HPDEALLC(PTworkl  ,IERR,1)
      write(nulout,*) 'IERR=',ierr
      CALL HPDEALLC(PTworkd  ,IERR,1)
      write(nulout,*) 'IERR=',ierr
      CALL HPDEALLC(PTd      ,IERR,1)
      write(nulout,*) 'IERR=',ierr
      CALL HPDEALLC(PTresid  ,IERR,1)
      write(nulout,*) 'IERR=',ierr
      CALL HPDEALLC(PTax     ,IERR,1)
      write(nulout,*) 'IERR=',ierr
      CALL HPDEALLC(PTselect ,IERR,1)
      write(nulout,*) 'IERR=',ierr
      CALL HPCHECK(IERR)
      IF(IERR.EQ.0)THEN
        WRITE(NULOUT,FMT='(5X,10("*")," No memory overflow"
     S        ," detected ",10("*"))')
      ELSE
        write(NULOUT,*) 'PROBLEM!!!!!'
      ENDIF
      write(nulout,*) 'after  dealloc'
      call vflush(nulout)
c
      end
c
c ------------------------------------------------------------------
c     matrix vector subroutine
c
c     Computes w <--- OP*v, where OP is the Hessian of Jo in space of
c     control vector
c ------------------------------------------------------------------

      SUBROUTINE AVSV (n, v, w) 3,33
      IMPLICIT NONE
#include "comdim.cdk"
#include "comlun.cdk"
#include "comcva.cdk"
#include "comgd0.cdk"
#include "comsv.cdk"
c
c Arguments
c
      integer n
      Real*8 v(n), w(n)
c
c Local variables
c
      integer i,JK,JLON,JLAT,numseg
      LOGICAL Ltalktome
      real*8 penergy1
      real*8, allocatable :: tlmstates(:,:,:,:)
c
      Ltalktome=.false.
      write(nulout,*) 'CSVNORM=',csvnorm,' NUMSEG=',numseg
c
      allocate(tlmstates(NI,NKGDIM,NJ,numseg))
c
      IF(CSVNORM.eq.'B') THEN
        IF(NPRECON.gt.0) THEN
          CALL pasqrt(v ,v ,n,nvadim)
        ENDIF
        CALL CAIN(NVADIM,V)
        CALL SPA2SP
        CALL SPGD
      ENDIF
      IF(CSVNORM.eq.'E') THEN
        CALL CAINGD(NSVDIM,V)
        CALL APPLYENERGY(.true.,.true.,.true.,.true.)
        CALL TRANSFER('GD10')
      ENDIF

      if(Ltalktome) then
      call calcenergy(penergy1)
      write(nulout,*) 'BEFORE TLM:',penergy1
      call vflush(nulout)
      endif

      nsim3d = nsim3d + 1
      write(nulout,*) 'CALLING TLM+ADJ:',nsim3d
      call vflush(nulout)
c
      do i=1,numseg
        write(nulout,*) 'CALLING TLM SEGMENT:',i
        call vflush(nulout)
        call putdx2('F')
        call getdx('F')
        do jlat=1,NJ
          do jk=1,NKGDIM
            do jlon=1,NI
              tlmstates(jlon,jk,jlat,i)=GD(jlon,jk,jlat)
            enddo
          enddo
        enddo
      enddo

      if(Ltalktome) then
      call calcenergy(penergy1)
      write(nulout,*) 'AFTER TLM:',penergy1
      call vflush(nulout)
      endif

c apply the final time norm to all saved states
      do i=1,numseg
        do jlat=1,NJ
          do jk=1,NKGDIM
            do jlon=1,NI
              GD(jlon,jk,jlat)=tlmstates(jlon,jk,jlat,i)
            enddo
          enddo
        enddo
        call applyenergy(.false.,.false.,.true.,.true.)
        call transfer('GD10')
c
        if(Ltalktome) then
        call calcenergy(penergy1)
        write(nulout,*) 'AFTER NORM:',penergy1
        call vflush(nulout)
        endif
c
        do jlat=1,NJ
          do jk=1,NKGDIM
            do jlon=1,NI
              tlmstates(jlon,jk,jlat,i)=GD(jlon,jk,jlat)
            enddo
          enddo
        enddo
      enddo

      call transfer('ZGD0')
      do i=numseg,1,-1
        write(nulout,*) 'CALLING ADJ SEGMENT:',i
        call vflush(nulout)
        do jlat=1,NJ
          do jk=1,NKGDIM
            do jlon=1,NI
              GD(jlon,jk,jlat)=GD(jlon,jk,jlat)+tlmstates(jlon,jk,jlat,i)
            enddo
          enddo
        enddo
        call putdx2('A')
        if(i.eq.1) then
          call getdx('A')
        else
          call getdx('X')
        endif
      enddo

      if(Ltalktome) then
      call calcenergy(penergy1)
      write(nulout,*) 'AFTER ADJ:',penergy1
      call vflush(nulout)
      endif

      IF(CSVNORM.eq.'E') THEN
        CALL APPLYENERGY(.true.,.true.,.true.,.true.)
        CALL TRANSFER('GD10')
        DO i=1,N
          W(i)=0.0
        ENDDO
        CALL CAINGDAD(NSVDIM,W)
      ENDIF
      IF(CSVNORM.eq.'B') THEN
        DO i=1,N
          VAZG(i) = 0.0
        ENDDO
        CALL SPGDA
        CALL SPA2SPAD
        CALL CAINAD(NVADIM,VAZG)
        IF(NPRECON.gt.0) THEN
          CALL pasqrt(vazg ,w ,nvadim,n)
        ELSE
          DO i=1,N
            W(i)=0.0
          ENDDO
          DO i=1,NVADIM
            W(i)=VAZG(i)
          ENDDO
        ENDIF
      ENDIF
c
      deallocate(tlmstates)
c
      return
      end