!-------------------------------------- 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 --------------------------------------
!
*DECK
SUBROUTINE VERMOD ( KULOUT ) 1,4
#if defined (DOC)
C
***s/r VERMOD - Computes the equivalent depths and vertical normal modes
C
C Author: Clive Temperton RPN/AES - Apr 1985
C
C Revision: Luc Fillion RPN/AES - MAY 93
C - Major modifications to suit the
C variational analysis environment.
C LAPACK was also introduced instead of NAG.
C
C Revision: Luc Fillion ARMA/AES - Jan 95
C - Minor modifications to clean the code.
C
C Revision: Luc Fillion - ARMA/AES - Oct 96
C - Strict use of VLEV and nflev and for CVCORD.ne.'PRESS'
* S. Pellerin *ARMA/SMC May 2000
* - Logical unit cleanup
* JM Belanger CMDA/SMC Aug 2000
* . 32 bits conversion
* (- Use of double precision BLAS routines DGEEV,
* DGETRI, DGETRF.
* - REAL*8 eps argument to MINV
* - MXMA8 instead of MXM
* - REAL*8 interface VIPSORT to IPSORT)
*
C----------
C
C Purpose: VERMOD computes the normal modes to be used by the variational analysis
C scheme and is compatible with the global spectral SEF model (P formulation)
C
C Arguments:
C i : KULOUT = Unit used for optional printing
C
C Documentation:
C DGEEV: LAPACK'S user's guide, SIAM 1992, p.14-15,136.
C DLINV is the matrix to be diagonalized.
C DLMAT contains the right-eigenvectors on output of DGEEV.
C DLDEP will contain the real part of the eigenvalues.
C N.B.: ON EXIT, DLINV HAS BEEN OVERWRITTEN.
C
C---------------------------------------------------------------------------
C
#endif
IMPLICIT NONE
C
*implicits
#include "pardim.cdk"
#include "comcst.cdk"
#include "comdim.cdk"
#include "comgem.cdk"
#include "comode.cdk"
#include "comlun.cdk"
#include "comlunla.cdk"
C
INTEGER KULOUT
C
INTEGER JI, JJ, ILEN, IERR, JLEV, JM, INFO, ISYM, IORD, ione
INTEGER INC, JG, ING, INR, JK, JP,iulmode,fnom,fclos
REAL*8 ZGAM0, ZDET, ZCON, ZOMAX, DUM, eps
INTEGER IPOINT(NFLEV)
REAL*8 ZWKV1(NFLEV), ZWKV2(NFLEV,3)
REAL*8 ZMAT1(NFLEV,NFLEV)
REAL*8 DLDEP(NFLEV), DLINV(NFLEV,NFLEV)
REAL*8 DLMAT(NFLEV,NFLEV)
REAL*8 DLW5(5*NFLEV), DLWM(NFLEV,NFLEV)
C
POINTER (PXIPOINT,IPOINT)
POINTER (PXWKV1,ZWKV1), (PXWKV2,ZWKV2)
POINTER (PXDLDEP,DLDEP), (PXDLINV,DLINV)
POINTER (PXDLMAT,DLMAT), (PXDLW5,DLW5)
POINTER (PXDLWM,DLWM), (PXZMAT1,ZMAT1)
C
*modules
EXTERNAL FEVOP, MINV, MXMA8, DGEEV, VIPSORT, DGETRF, DGETRI
EXTERNAL HPALLOC, HPDEALLC, ABORT3D
*
**
C
C*0 . ALLOCATE LACAL ARRAYS
C ---------------------
C
ILEN = NFLEV
CALL HPALLOC ( PXIPOINT,MAX(ILEN,1),IERR,8)
CALL HPALLOC ( PXWKV1,MAX(ILEN,1),IERR,8)
CALL HPALLOC ( PXDLDEP, MAX(ILEN,1), IERR, 8 )
C
ILEN = 3*NFLEV
CALL HPALLOC ( PXWKV2,MAX(ILEN,1),IERR,8)
C
ILEN = 5*NFLEV
CALL HPALLOC ( PXDLW5, MAX(ILEN,1), IERR, 8 )
C
ILEN = NFLEV*NFLEV
CALL HPALLOC ( PXZMAT1, MAX(ILEN,1), IERR, 8 )
CALL HPALLOC ( PXDLMAT, MAX(ILEN,1), IERR, 8 )
CALL HPALLOC ( PXDLINV, MAX(ILEN,1), IERR, 8 )
CALL HPALLOC ( PXDLWM, MAX(ILEN,1), IERR, 8 )
C
C
C*1 VERTICAL PROBLEM: COMPUTATION OF THE EQUIVALENT DEPTHS AND VERTICAL MODES
C -------------------------------------------------------------------------
C
eps=0.0D0
100 CONTINUE
IF(NFLEV.GT.1) THEN
WRITE(KULOUT,1000)
ZGAM0 = RKAPPA * TSTAR
C
DO 110 JLEV = 1 , NFLEV-1
ZWKV1(JLEV) = VLEV(JLEV+1) - VLEV(JLEV)
110 CONTINUE
ZWKV1(NFLEV) = 0.0
CALL FEVOP
( ZMAT1, VMAT, ZWKV1, ZWKV2, ZGAM0, NFLEV )
CALL MINV
( ZMAT1, NFLEV, NFLEV, VINV, ZDET, eps, 0, +1 )
CALL MXMA8( ZMAT1, 1, NFLEV, VMAT, 1, NFLEV, VINV, 1, NFLEV,
$ NFLEV, NFLEV, NFLEV)
C
C TRANSFER THE RESULT TO DOUBLE PRECISION PRIOR TO LAPACK SOLVER
C
DO 130 JJ = 1, NFLEV
DO 120 JI = 1, NFLEV
DLINV(JI,JJ) = VINV(JI,JJ)
120 CONTINUE
130 CONTINUE
C
C COMPUTE THE EIGENVALUES/EIGENVECTORS
C
CALL DGEEV ( 'N', 'V', NFLEV, DLINV, NFLEV, DLDEP,
+ DLW5, DUM, 1, DLMAT, NFLEV, DLW5(NFLEV+1),
+ 4*NFLEV, INFO )
C
C THE COMPUTED EIGENVECTORS HAVE BEEN NORMALIZED BY DGEEV
C TO HAVE EUCLIDEAN NORM EQUAL TO 1 AND LARGEST COMPONENT REAL.
C MAKE SURE ALL EIGENVALUES ARE REAL
C
DO 140 JLEV = 1, NFLEV
IF(DLW5(JLEV).NE.0.) THEN
WRITE(KULOUT,*)'*********************************************'
WRITE(KULOUT,*)'IMAG PART EIGNVALUE FOR VRTCL PRBLM NOT ZERO '
WRITE(KULOUT,*)'PROGRAM STOPS'
WRITE(KULOUT,*)'********************************************'
CALL ABORT3D
(NULOUT,'ABORT IN VERMOD')
ENDIF
140 CONTINUE
C
C SORT THE EIGENVALUES
C
DO 150 JLEV = 1, NFLEV
150 ZWKV1(JLEV) = DLDEP(JLEV)
CALL VIPSORT
( IPOINT, ZWKV1, NFLEV )
C
C DO A PERMUTATION OF THE EIGENVECTORS/VALUES
C
DO 170 JLEV = 1, NFLEV
JP = IPOINT(JLEV)
DO 160 JK = 1, NFLEV
160 DLWM(JK,JLEV) = DLMAT(JK,JP)
EQDEPTH(JLEV) = DLDEP(JP)
170 CONTINUE
C
C SET VMAT AND PREPARE FOR INVERSION
C
DO 190 JLEV = 1, NFLEV
DO 180 JK = 1, NFLEV
VMAT(JK,JLEV) = DLWM(JK,JLEV)
DLINV(JK,JLEV) = DLWM(JK,JLEV)
180 CONTINUE
190 CONTINUE
C
C COMPUTE THE INVERSE OF DLMAT: RESULT IN DLINV
C N.B. DLINV MUST CONTAIN DLMAT
C
CALL DGETRF ( NFLEV, NFLEV, DLINV, NFLEV, IPOINT, INFO )
CALL DGETRI ( NFLEV, DLINV, NFLEV, IPOINT, DLWM,
+ NFLEV*NFLEV, INFO )
C
C CONVERT EIGENVALUES INTO EQUIVALENT DEPTHS
C
ZCON = RD / RG
DO 195 JLEV = 1 , NFLEV
EQDEPTH(JLEV) = ZCON / EQDEPTH(JLEV)
195 CONTINUE
ELSE
EQDEPTH(1) = BRTPHGT
ENDIF
C
C*2 WRITE EQUIVALENT DEPTHS AND VERTICAL MODES ON FILE
C --------------------------------------------------
C
200 CONTINUE
C MOVE NORMAL MODE RESULTS FROM REAL*8 TO REAL TO BE USED
C BY THE VARIATIONAL ANALYSIS SCHEME OR SEF MODEL.
C
DO 220 JJ = 1, NFLEV
DO 210 JI = 1, NFLEV
C N.B. VMAT WAS PREVIOUSLY SET
VINV(JI,JJ) = DLINV(JI,JJ)
210 CONTINUE
220 CONTINUE
C
IF(LWRTMOD) THEN
iULMODE = 0
close(nutemp)
open (unit=nutemp,file='vertical_modes.od')
ione = 1
write(nutemp,910) nflev,nflev,ione
write(nutemp,920) ione
do ji = 1,nflev
write(nutemp,900) (vmat(ji,jj),jj=1,nflev)
enddo
close(nutemp)
900 format(100(E13.7,3X))
910 format(3(I4,1X))
920 format(I4)
ENDIF
C
C PRINT EQUIVALENT DEPTHS ON STANDARD OUTPUT
C
DO 230 JLEV = 1 , NFLEV
WRITE(KULOUT,1040) JLEV, EQDEPTH(JLEV)
230 CONTINUE
C
C*3 DEALLOCATE LOCAL ARRAYS
C -----------------------
C
300 CONTINUE
CALL HPDEALLC ( PXIPOINT, IERR, 1 )
IF(IERR.NE.0)WRITE(NULOUT,*)'OVERFLOW IN IPOINT, IERR =',IERR
CALL HPDEALLC ( PXWKV1, IERR, 1 )
IF(IERR.NE.0)WRITE(NULOUT,*)'OVERFLOW IN ZWKV1, IERR =',IERR
CALL HPDEALLC ( PXDLDEP, IERR, 1 )
IF(IERR.NE.0)WRITE(NULOUT,*)'OVERFLOW IN DLDEP, IERR =',IERR
CALL HPDEALLC ( PXWKV2, IERR, 1 )
IF(IERR.NE.0)WRITE(NULOUT,*)'OVERFLOW IN ZWKV2, IERR =',IERR
CALL HPDEALLC ( PXDLW5, IERR, 1 )
IF(IERR.NE.0)WRITE(NULOUT,*)'OVERFLOW IN DLW5, IERR =',IERR
CALL HPDEALLC ( PXZMAT1, IERR, 1 )
IF(IERR.NE.0)WRITE(NULOUT,*)'OVERFLOW IN ZMAT1, IERR =',IERR
CALL HPDEALLC ( PXDLMAT, IERR, 1 )
IF(IERR.NE.0)WRITE(NULOUT,*)'OVERFLOW IN DLMAT, IERR =',IERR
CALL HPDEALLC ( PXDLINV, IERR, 1 )
IF(IERR.NE.0)WRITE(NULOUT,*)'OVERFLOW IN DLINV, IERR =',IERR
CALL HPDEALLC ( PXDLWM, IERR, 1 )
IF(IERR.NE.0)WRITE(NULOUT,*)'OVERFLOW IN DLWM, IERR =',IERR
C
1000 FORMAT(/2X,' COMPUTATION OF THE VERTICAL NORMAL MODES'/)
1020 FORMAT(/' PROBLEM, BAD ORDER MATRICES IN VERMOD')
1030 FORMAT(/' PROBLEM WITH LAPACK SUBROUTINE , IERR=',
1 I3,' ...STOP...')
1040 FORMAT(/10X,' MODE',I3,' EQUIVALENT DEPTH=',F12.4)
C
C
RETURN
END