exci/res_exci_bcm.f

Fortran project RESTRAX, source module exci/res_exci_bcm.f.

Source module last modified on Thu, 11 May 2006, 18:33;
HTML image of Fortran source automatically generated by for2html on Mon, 29 May 2006, 15:06.


#////////////////////////////////////////////////////////////////////////////////
#////
#////  R E S T R A X   4.8.1    EXCI
#////
#//// Subroutine called by RESTRAX to get values of excitation energy (OMEXC)
#//// and scattering cross-section (SQOM) for given QHKL,E values stored in Q(i)
#//// Permits to define up to 6 different branches of S(Q,E)
#////
#//// You can use this file as a template
#//// Refer to DON'T CHANGE .. END blocks for the code to be preserved
#////
#//// J. Saroun (saroun@ujf.cas.cz) , March 2005
#//// Read attached documentation or visit RESTRAX home page for help: 
#//// http://omega.ujf.cas.cz/restrax
#////////////////////////////////////////////////////////////////////////////////
#
#                           ***  ARGUMENTS ***
# input:
# Q(1:4)   ... (H,K,L,E) values 
# ICOM<-10 ... initialization (called only once when loaded at runtime )
# ICOM=0   ... initialization (run usually before each  [M]FIT or INIT  commands)
# ICOM=-1  ... only excitation energies are used (e.g. for plotting disp. branches)
# ICOM=-2  ... only S(Q,E) values are used (e.g. for plotting S(Q,E) maps)
# ICOM>0   ... should return both excitation energies and S(Q,E). ICOM=index of supplied event
#
# output:
# OMEXC(1:6) ... excitation energies for 1..nbr branches for Qhkl = Q(1:3) 
# SQOM(1:6)  ... S(Q,E) values for 1..nbr branches
#
#                           *** SHARED DATA ***
#
# Following fileds are available via common variables declared in the *.inc files:
#
# Monte Carlo ray-tracing results:
#-----------------------------------
# accessible only if ICOM>=0 !!
# REAL*4 QOM(1:4,j),PQOM(j) .... value of (Q,E) and weight for j-th event
# IQOM(j) .... the index of data set corresponding to given j-th event.
# NQOM(k) .... partitioning of the QOM, PQOM ... arrays, i.e. the number
#              of events stored for the k-th data set is NQOM(k)-NQOM(k-1).
# NDATQOM .... index of actual data set, for which the scan profile is accumulated
#              Use this index to define specific free parameters for different data sets 
#
# Instrument setting:
#-----------------------------------
# REAL*4 QOM0(1:4,k)   .... Spectrometer position (Q,E) for k-th data set
#
# Unit vectors in rec. lat. units:
#-----------------------------------
# REAL*8 PARAM(1:MPAR)         ... free model parameters
# INETEGR*4 FIXPARAM(1:MPAR)   ... fixed parameters. Set FIXPARAM(i)=0 to make  
#                                  the i-th parameter fixed)
# NTERM                        ... number of free model parameters (<=64)
# NBR                          ... number of branches defined by EXCI (<=6)
# REAL*8 WEN(1:6)              ... widths of the disp. branches. 
# CHARACTER*10 PARNAME(1:MPAR) ... names of free parameters
#                         
# Outside EXCI, WEN is used only as a flag to check, whether scattering 
# is difuse (WEN>0) or not (WEN=0). The convolution method is selected according 
# to this flag.  
#
#                           *** SHARED SUBROUTINES ***
#                        (see source files for details)
#   in this module:
#   SUBROUTINE READEXCIPAR     ... Read initial values of model variables 
#   in exci_io.f:
#   SUBROUTINE SETEXCIDEFAULT  ... Set default values to common EXCI variables 
#   in reclat.f:
#   SUBROUTINE POLVECT(Q,TAU,SIG1,SIG2,SIG3,ICOM) ... Get polarization unit vectors with 
#                                                     respect to q=TAU-Q 
#   REAL*8 FUNCTION QxQ(A,B)   ... Scalar product of vectors A,B in non-carthesian rec. lattice coordinates
#   SUBROUTINE QNORM(X,QRLU,QANG)   ... Norm of a vector X in non-carthesian rec. lattice coordinates
#///////////////////////////////////////////////////////////////////////////////////
#------------------------------------------------------------------------------
      SUBROUTINE EXCI(icom,q,omexc,sqom)
# Bond charge model: phonons in diamond lattice (Si, Ge, ...)
#     REFERENCE: W.WEBER, PHYS.REV.B,VOL.15,NO.10,(1977),4789.
#                O.H.NIELSEN, W.WEBER, COMPUTER PHYSICS COMM.(1979)
# The energy profile S(E) is either zero-width [wen(j)=0], or Lorenzian [wen(j)>0]
#------------------------------------------------------------------------------
      implicit none
      
#----------------------- *** DON'T CHANGE *** ------------------------------
      INCLUDE 'const.inc'
#      INCLUDE 'inout.inc'
      INCLUDE 'exci.inc'
      integer*4 icom,excinit
      real*8 q(4),omexc(6),sqom(6)
      
#-------------------------- *** END *** ------------------------------------
      
# **** Local user declarations ****

      real*8 pihalf,tau(3)
      parameter (pihalf=1.570796327d0)   
      character*16 cname(6)
      integer*4 j,i,k,nptmc,nlog
      real*8 bf,kt,z,arg
      complex eivec(mqom,6,6),eivbmc(6,6),f(6),phfac
      real*4  qphon(3),omegaint(mqom,6),ombmc(6) !,phase
      real*8 lastchk
      integer*4 lastimod
      save eivec,omegaint
      data tau/2.d0,2.d0,0.d0/
#// internal variables describing the model
#// This common is not shared with the rest of RESTRAX, only within this file. 
      integer*4 imodel
      real*8 temp,fqsq(6),omexc0(6),omscal
      common /excipar/ cname,temp,fqsq,omexc0,omscal,imodel
#// for debugging only, use with bcm_debug.f:      
#      real*4 q2(3),om2(6)
#      real*8 f2(6)
#      complex cvec(6),eiv2(6,6)
      
      
#------------------- *** DATA section *** ----------------

# **** DEFAULT values of internal model variables ****
      data lastchk/1.1d0/ 
      data lastimod/-1/
      data cname / 'Diamond', 'Silicon', 'Germanium', 'alpha-Tin', 'Ge70',
     &            'Ge76'/     

      data temp/300.d0/
      data fqsq/6*1.d0/
#      DATA tau/1.,0.,0./tauread/0/
      data omscal/1.0/
      data omexc0/5*0.0,5.0/
      data excinit/0/

#***********************************************************************************
# MODEL INITIALIZATION (ICOM<-10)
#***********************************************************************************
#-- called only once when loaded at runtime 
#-- set some values shared with RESTRAX if different from default

      if (icom.lt.-10) then

# Set model identification string: 
        phontitle=
      'Bond charge model: phonons in diamond lattice (Si, Ge, ...)' 
              
#// Define fixed parameters (=0), default: all free (=1)
#        FIXPARAM(1)=0  ! let Intensity fixed !!

#// Number of branches ****    
        nbr=6

#// Initial widths in energy, default=1meV 
#// Set wen(i)=0 for zero-width branches
        do i=1,6
          wen(1)=0.d0 
        enddo  

#**** How to read file with parameters (default=1):
#**** (0) never (1) at program start or on INIT command (2) each time MFIT is called   
#        EXCREAD=0 
      
# Set name of file with model parameters (if different from default exc.par)
        phonname= 'bcm.par'
      
# Define names of free parameters for i>2:
        parname(3)= 'EN_scale'
        parname(4)= 'W_LO'       
        parname(5)= 'W_TO1'       
        parname(6)= 'W_TO2'       
        parname(7)= 'W_LA'       
        parname(8)= 'W_TA1'       
        parname(9)= 'W_TA2'          

#// number of free model parameters
        nterm=9
        
        write(*,*)  'EXCI: set default'
        return
      endif
 
#----------------------- *** DON'T CHANGE *** ------------------------------  
      if ((icom.ne.0).and.(excinit.ne.0)) goto 1      
#---------------------------- *** END *** ----------------------------------

#***********************************************************************************
# MODEL INITIALIZATION (ICOM=0)                    
#***********************************************************************************
#-- called before each [M]FIT or INIT command  
#

# *** initialize BCM
      call parbcm(imodel,0)      

      nptmc=nqom(mdat)  ! get number of events stored for all data sets
      
# *** tabulate BCM if model or QOM values changed
      if (lastimod.ne.imodel.or.lastchk.ne.chkqom) then
10    format( 'Tabulating BCM frequencies & eigenvectors for ',i8, ' points',$)
11    format( '.',$)
      write(*,10) nptmc
      nlog=nptmc/20
      do i=1,nptmc
        if (nlog.gt.0.and.mod(i,nlog).eq.0) write(*,11)
        do j=1,3
          qphon(j) = qom(j,i)  !-tau(j)
        enddo
        call bcmres(qphon,ombmc,eivbmc)
        do k=1,6
          omegaint(i,k) = ombmc(k)
          do j=1,6
            eivec(i,j,k) = eivbmc(j,k)  !k ... branch, j ... comp.
          enddo
        enddo
      enddo
      write(*,*) 'finished.'
      lastchk=chkqom
      lastimod=imodel
      endif

#// Assign free model parameters to param():
      param(3)=omscal
      do i=1,nbr 
        param(3+i)= wen(i) 
      enddo  
      
#----------------------- *** DON'T CHANGE *** -------------------------

      excinit=1
      return
#---------------------------- *** END *** ----------------------------------
1     continue      

#********************************************************************************
#                                                                   
#                   EXECUTION PART (ICOM<>0)                        
#
# This part is called many times during the fitting procedure 
# =>  should be as fast as possible
# 
#// Do whatever you want in the following code. 
#// EXCI MUST RETURN: 
#// OMEXC(i) ... excitation energies for first NBR branches (i=1..6)
#// SQOM(i)  ... dS/dOmega/dE 


# ICOM=-1 => only OMEXC(i) values are used by RESTRAX to plot the branches.
# ICOM=-2 => both OMEXC and SQOM are needed, but no data set is provided (used for mapping S(Q,E)
# Otherwise, ICOM refers to the event number in the QOM array
#   => ICOM can be used e.g. as an index to internal lookup tables of EXCI etc...                                
#********************************************************************************
      

#------------------!! OBLIGATORY !!-------------------------

#// Assign values in the PARAM array to the local model variables
#// if you don't work with the PARAM() array directly
#// REMEMBER: PARAM(1,2) are reserved for Scale and Background

      omscal = param(3)
      do j=1,nbr
        wen(j) = param(3+j)
      enddo

#// don't allow wen->0, if it is not a zero-width branch  !!      
      do i=1,nbr
         wen(i)=abs(wen(i))
         if(wen(i).ne.0.d0.and.wen(i).lt.1e-3) wen(i) = 1.e-3 
      enddo     
        
#----------------------!! END !! -------------------------

#*  If ICOM=-1, return only excitation energies
      if(icom.eq.-1) then
        do j=1,3
          qphon(j) = q(j)  !-tau(j)
        enddo
        call bcmres(qphon,ombmc,eivbmc)
        do j=1,6
          omexc(j) = omscal*ombmc(j)
          sqom(j) = 1.
       enddo
       return
      endif

#///  Returns up to six energies and cross-sections:
         
      do i=1,3
         qphon(i) = q(i)  !-tau(i)  ! get actual phonon q
      enddo
      
#/// get energies and eigenvectors from table if possible
#/// otherwise call bcmres

#// check the difference between tabulated and actual Q-value      
      z=0.d0
      if (icom.gt.0) then
        do j=1,3
         z=z+abs(qom(j,icom)-q(j))
        enddo 
      endif
        
#// calculate eigenvalues/vectors if qom differes from initialization input or if icom<0
      if (icom.le.0.or.z.gt.1e-5) then
        call bcmres(qphon,ombmc,eivbmc)      
      else
#// otherwise use values calculated during initialization
        do k=1,nbr
          ombmc(k) = omegaint(icom,k)
          do j=1,6
            eivbmc(j,k)=eivec(icom,j,k)   !k ... branch, j ... comp.
          enddo
        enddo
      endif   

#// scale energies               
      do k=1,nbr
        omexc(k) = omscal*ombmc(k)
      enddo
 
# *** calculate cross-section

#      phase = pihalf*(tau(1)+tau(2)+tau(3))
#      phfac = cexp(cmplx(0.D0,phase))
      phfac=cmplx(1.0,0.0)
      do k=1,nbr
        f(k) = cmplx(0,0)
          do j=1,3
            f(k)=f(k)+q(j)*(eivbmc(j,k)*phfac+eivbmc(j+3,k))
          enddo
        fqsq(k)=(abs(f(k)))**2
      enddo


# This was for debugging  only, comment the following include out for normal use
#-----------------------------------------------------------------------
#      INCLUDE 'bcm_debug.f'  
#-----------------------------------------------------------------------
      
           
100   kt=temp/11.609  ! 11.609 ... conversion kT -> meV
      do k=1,nbr
#// bf = Bose factor x omega
        z=exp(-omexc(k)/kt)  
        if (abs(1.0-z).lt.1.d-5) then
          bf=kt
        else
          bf=omexc(k)/(1-z)
        endif          
        if(fqsq(k).eq.0) then
          sqom(k)=0.
        else
          if(abs(omexc(k)).lt.0.0001) omexc(k)=0.0001
          if(wen(k).gt.0.) then  ! apply finite width
#             arg = 2.*(omexc(k)-q(4))/wen(k)
#             sqom(k) = bf*fqsq(k)*2./(wen(k)*(1.+arg**2))
# modified: use damped oscillators instead of Lorentz dist.          
            arg=((q(4)**2-omexc(k)**2)**2+wen(k)**2*q(4)**2)/wen(k)
            if(abs(arg).lt.1.e-8) arg=1.e-8
            sqom(k) = bf*fqsq(k)/abs(arg)
           else
             if(abs(omexc(k)).lt.1.e-4) omexc(k)=1.e-4
             sqom(k) = bf*fqsq(k)/omexc(k)**2
           endif
        endif
      enddo

      end
      
#------------------------------------------------------------------------------
      SUBROUTINE REPEXCIPAR
# REPORT model ID and input parameters as needed
#------------------------------------------------------------------------------
      implicit none
      INCLUDE 'const.inc'
#      INCLUDE 'inout.inc'
      INCLUDE 'exci.inc'
      
      character*16 cname(6)
      integer*4 imodel
      real*8 temp,fqsq(6),omexc0(6),omscal
      common /excipar/ cname,temp,fqsq,omexc0,omscal,imodel
      
      write(*,*)  'EXCI: ',trim(phontitle)
#// Report some model values:  
      write(*,*)  'Material: ',trim(cname(imodel))
12    format( ' Temperature [K]: ',g10.4)
      write(*,12) temp
     
      end
      
#------------------------------------------------------------------------------
      SUBROUTINE READEXCIPAR
# Read values of model variables used by EXCI
# Call by RESTRAX when requiared
# File is opened and closed by RESTRAX, don't call OPEN/CLOSE here !!!
#------------------------------------------------------------------------------
      implicit none
      INCLUDE 'const.inc'
#      INCLUDE 'inout.inc'
      INCLUDE 'exci.inc'
     

      character*16 cname(6)
      integer*4 imodel
      real*8 temp,fqsq(6),omexc0(6),omscal
      common /excipar/ cname,temp,fqsq,omexc0,omscal,imodel
      integer*4 i
           
      rewind(excunit)  ! call rewind for compatibility with g77
# read some model parameters from file 
      read (excunit,*,err=998) temp   ! temperature     
      if (temp.le.0.) temp=0.01

      read (excunit,*,err=998) imodel   ! Diamond, Si, Ge, ...
      if (imodel.gt.6.or.imodel.le.0) imodel=1
      read (excunit,*,err=30) (wen(i),i=1,6)  ! # widths (optional)              
            
30    continue
      write(*,*)  'Parameters updated from '//phonname
      return
      
998   write(*,*)  'Format error?! Cannot read excitation parameters.'
      return
      end
      
      
#--------------------------------------------------------------------
      SUBROUTINE BCMRES (q,frmev,vec)
#     Bond charge model - modified by J.K. Oct-1998
#
# ***  procedure for Restrax calls
# *** PARAM is called PARBCM to avoid confusion
#--------------------------------------------------------------------

      real q(3),frmev(6),wr(6,6),wi(6,6)
      complex vec(6,6)
#
#     BOND CHARGE MODEL IN DIAMOND-SYMMETRY
#     REFERENCE: W.WEBER, PHYS.REV.B,VOL.15,NO.10,(1977),4789.
#                O.H.NIELSEN, W.WEBER, COMPUTER PHYSICS COMM.(1979)
#
#     INPUT:
#
#     Q ... WAVEVECTOR
#           WE DEFINE THE 1BZ BY:
#           X=(1,0,0)
#           L=(.5,.5,.5)
#
#     OUTPUT:
#
#     FRTHZ ...  ARRAY CONTAINING THE EIGENVALUES IN INCREASING ORDER
#     WR,WI ...  REAL & IMAGINARY PART OF THE CORRESPONDING EIGENVECTORS
#                FOR FRTHZ(I) THE EIGENVECTOR IS:
#                W(J,I),J=1..6
#
#                COMMON /PRAM/  from param.fortran
#
#     FCBC ..... BOND CHARGE FORCE CONSTANTS, USING NOTATION OF W.WEBER:
#                1-ALFA; 2-BETA; 3-MU; 4-NU; 5-LAMBDA; 6-DELTA
#     FCION .... ION-ION FORCE CONSTANTS, SAME AS FCBC, BUT WITH A PRIME
#                IN W.WEBER'S NOTATION.
#     MASS ..... MASS OF THE ATOM.
#     SCALE .... A SCALE FACTOR TO CONVERT FREQUENCIES INTO THZ, SINCE
#                INTERNALLY WE USE FORCE CONSTANT UNITS E**2/VA.
#     ACELL2 ... 1/2 OF CELL LENGTH, DENOTED R0 IN W.WEBER'S ARTICLE.
#     Z2EPS .... Z**2/EPSILON (MODEL PARAMETER).
#
#     COMMON AREA OUTPUT, WHICH THE USER MIGHT WISH TO USE:
#
#     EXPRAM ... EXPERIMENTAL RAMAN FREQUENCY.
#     THRAM .... THEORETICAL RAMAN FREQUENCY.
#     CELAST ... ELASTIC CONSTANTS ARRANGED AS:
#                C11  C12  C44  EXPERIMENTAL
#                C11  C12  C44  CALCULATED
#
      common /pram/ fcbc(6),fcion(6),mass,scale,acell2,z2eps,expram,thr
     1am,celast(3,2)
      real mass
#
#     QPI ... WAVEVECTOR MULTIPLIED BY PI
#
      common /qvect/ qpi(3)
#
#     A ... PRIMITIVE TRANSLATIONAL VECTORS.
#           A(I,J) IS COORDINATE I OF THE J'TH VECTOR.
#     R12 . CONTAINS ALL THE RELATIVE COORDINATES OF THE IONS AND
#           BOND CHARGES, I.E. X(K')-X(K).
#           R12(I,J) IS COORDINATE J OF RELATIVE VECTOR I.
#
      common /basis/ a(3,3),r12(12,3)
#
      real crr(6,6),cri(6,6),ctr(6,12),cti(6,12),cs(12,12),tr(6,12),ti(6
     1,12),s(12,12),ur(12,6),ui(12,6),tst1(6,6),tst2(6,6),tstr(6,6),tsti
     2(6,6),rr(6,6),ri(6,6)
              
     
      data pi/3.1415926535898/
      data clight,mev/33.3564095198,4.135701/
#
#     AT GAMMA, WE CHOOSE Q TO LIE ALONG (100)
      if (abs(q(1))+abs(q(2))+abs(q(3)).le.1.0e-6) q(1)=0.00001
#
      do i=1,3
        qpi(i)=q(i)*pi       !     WE NEED THE WAVEVECTOR TIMES PI
      end do
#
#     CONSTRUCT THE COULOMB MATRICES CR,CT,CS
      call CCMAT (crr,cri,ctr,cti,cs)
#
      call MATOP ( '.    ',crr,6,6,crr,6,6,crr,z2eps)
      call MATOP ( '.    ',cri,6,6,cri,6,6,cri,z2eps)
      call MATOP ( '.    ',ctr,6,12,ctr,6,12,ctr,z2eps)
      call MATOP ( '.    ',cti,6,12,cti,6,12,cti,z2eps)
      call MATOP ( '.    ',cs,12,12,cs,12,12,cs,z2eps)
#
#     CONSTRUCT MATRIX T
      call TMAT (tr,ti)
      call MATOP ( '+    ',tr,6,12,ctr,6,12,tr,1.0)
      call MATOP ( '+    ',ti,6,12,cti,6,12,ti,1.0)
      call MATOP ( 'TRANS',tr,6,12,ur,12,6,tr,1.0)
      call MATOP ( 'TRANS',ti,6,12,ui,12,6,ti,-1.0)
#
#     CONSTRUCT MATRIX S
      call SMATX (s)
#
      call MATOP ( '+    ',s,12,12,cs,12,12,s,1.0)
      call MATOP ( 'INV  ',s,12,12,cs,12,12,cs,1.0)
      call MATOP ( '*    ',tr,6,12,cs,12,12,ctr,1.0)
      call MATOP ( '*    ',ctr,6,12,ur,12,6,tst1,1.0)
      call MATOP ( '*    ',ti,6,12,cs,12,12,cti,1.0)
      call MATOP ( '*    ',cti,6,12,ui,12,6,tst2,1.0)
      call MATOP ( '-    ',tst1,6,6,tst2,6,6,tstr,1.0)
      call MATOP ( '*    ',cti,6,12,ur,12,6,tst1,1.0)
      call MATOP ( '*    ',ctr,6,12,ui,12,6,tst2,1.0)
      call MATOP ( '+    ',tst1,6,6,tst2,6,6,tsti,1.0)
#
#     CONSTRUCT MATRIX R
      call RMAT (rr,ri)
#
      call MATOP ( '+    ',rr,6,6,crr,6,6,rr,1.0)
      call MATOP ( '+    ',ri,6,6,cri,6,6,ri,1.0)
      call MATOP ( '-    ',rr,6,6,tstr,6,6,rr,1.0)
      call MATOP ( '-    ',ri,6,6,tsti,6,6,ri,1.0)
#
#     THE TOTAL DYNAMICAL MATRIX: (WR,WI)
#
#     WEIGHT WITH MASS FACTORS
      call MATOP ( '.    ',rr,6,6,wr,6,6,wr,1.0/mass)
      call MATOP ( '.    ',ri,6,6,wi,6,6,wi,1.0/mass)

# check eigenvalues
#10    format(a,6(2x,G14.8))
#11    format(a,6(2x,G10.4))
#      write(*,10) 'Q: ',(Q(I),I=1,3)
#      write(*,10) 'WI:',(WI(I,I),I=1,6)
#      do i=1,6
#      do j=1,6
#         crr(i,j)=WR(i,j)
#
#      enddo
#      enddo

#
#     SET THE (HOPEFULLY) SMALL DIAGONAL IMAGINARY ELEMENTS = 0
#

      do 120 i=1,6
      wi(i,i)=0.0
120   continue
#
#     SOLVE THE EIGENVALUE PROBLEM
#

      call MATOP ( 'EIG  ',wr,6,6,wi,6,6,tst1,1.0)

      
#      write(*,10) 'TESTING EISPACK: ',ierr*1.D0

#      do j=1,6
#        sum=0.
#
#
#
        
#        write(*,11) 'OM, NORM(e): ',tst1(j,j),sqrt(sum)
#      do i=1,6
#        tstr(i,j)=0.
#
#
#
#
#
#
#     &         wr(i,j),wi(i,j)
#      enddo
#      enddo
#      pause
      
      
#
#     TAKE "SQUARE ROOT"
#
      do 130 i=1,6
      frmev(i)=sqrt(abs(tst1(i,i)*scale))*sign(1.0,tst1(i,i))*mev
130   continue
#
# *** the eigenvectors are in C-convention
#
# *** to convert to Strauch's-convention
# *** (include trivial phase factor) uncomment the following ***
#
#     phase = .5*(qpi(1)+qpi(2)+qpi(3))
#      cph = cos(phase)
#      sph = sin(phase)
#      do j=1,6
#        do i=4,6
#          wwr =wr(i,j)*cph-wi(i,j)*sph
#          wwi =wi(i,j)*cph+wr(i,j)*sph
#          wr(i,j) = wwr
#          wi(i,j) = wwi
#        enddo
#      enddo
# *** down to here *********************************************
#
      do i=1,6            ! added by STZO
        do j=1,6
          vec(i,j)=cmplx(wr(i,j),wi(i,j))
        end do
      end do
#
      return
      end

#
#--------------------------------------------------------------------
      SUBROUTINE CCMAT (crr,cri,ctr,cti,cs)
      real crr(6,6),cri(6,6),ctr(6,12),cti(6,12),cs(12,12)
#
#     CONSTRUCT THE COULOMB MATRICES BY CALLING COULCF
#
#     OUTPUT:
#     CRR,CRI ... REAL & IMAG PART OF MATRIX CR
#     CTR,CTI ... REAL & IMAG PART OF MATRIX CT
#     CS ........ MATRIX CS (A REAL MATRIX BY SYMMETRY)
#
#--------------------------------------------------------------------
      complex bbb(3,3)
      real aa(6,6),ab(6,6),ba(6,6),bb(6,6),r(3),bbbr(12,3,3),bbbi(12,3,3
     1),a1(3,3),a2(3,3),a3(3,3),a4(3,3)
      common /qvect/ qpi(3)
      common /basis/ a(3,3),r12(12,3)
      data eps/1.3/,icount/0/,tp/83.99769/
      save a1,a2,a3,a4  ! static variables, added by J.S. 
#
#     CALL GENERAL COULOMB COEFFICIENT ROUTINE
#
      do 110 i=1,12
      do 100 j=1,3
      r(j)=r12(i,j)
100   continue
      call COULCF (bbb,qpi,a,r,eps)
#
      do 110 j=1,3
      do 110 k=1,3
      bbbr(i,j,k)=real(bbb(j,k))
      bbbi(i,j,k)=aimag(bbb(j,k))
110   continue
911   format(12f6.2)
#
#     GENERATE MATRIX CR
#
      do 120 i=1,3
      do 120 j=1,3
      i3=i+3
      j3=j+3
      crr(i,j)=4.0*bbbr(1,i,j)
      crr(i3,j3)=crr(i,j)
      cri(i,j)=4.0*bbbi(1,i,j)
      cri(i3,j3)=cri(i,j)
      crr(i,j3)=4.0*bbbr(2,i,j)
      crr(i3,j)=crr(i,j3)
      cri(i,j3)=4.0*bbbi(2,i,j)
      cri(i3,j)=-cri(i,j3)
120   continue
#
#     GENERATE MATRIX CT
#
      do 140 i=1,3
      do 130 j=1,3
      ctr(i,j)=-2.0*bbbr(3,i,j)
      cti(i,j)=-2.0*bbbi(3,i,j)
      ctr(i,j+3)=-2.0*bbbr(4,i,j)
      cti(i,j+3)=-2.0*bbbi(4,i,j)
      ctr(i,j+6)=-2.0*bbbr(5,i,j)
      cti(i,j+6)=-2.0*bbbi(5,i,j)
      ctr(i,j+9)=-2.0*bbbr(6,i,j)
      cti(i,j+9)=-2.0*bbbi(6,i,j)
130   continue
      do 140 j=1,12
      ctr(i+3,j)=ctr(i,j)
      cti(i+3,j)=-cti(i,j)
140   continue
#
#     GENERATE MATRIX CS
#
#     FIRST TIME CALLING CCMAT, WE BUILD THE MATRICES A1..A4
      if (icount.eq.1) go to 150
      icount=1
      call MAT33 (0.,-1.,-1.,-1.,0.,-1.,-1.,-1.,0.,tp,a1)
      call MAT33 (0.,1.,-1.,1.,0.,1.,-1.,1.,0.,tp,a2)
      call MAT33 (0.,1.,1.,1.,0.,-1.,1.,-1.,0.,tp,a3)
      call MAT33 (0.,-1.,1.,-1.,0.,1.,1.,1.,0.,tp,a4)
150   do 160 i=1,3
      do 160 j=1,3
      i3=i+3
      j3=j+3
      aa(i,j)=bbbr(1,i,j)+a1(i,j)
      aa(i3,j3)=bbbr(1,i,j)+a2(i,j)
      aa(i,j3)=bbbr(7,i,j)
      aa(i3,j)=aa(i,j3)
      bb(i,j)=bbbr(1,i,j)+a3(i,j)
      bb(i3,j3)=bbbr(1,i,j)+a4(i,j)
      bb(i,j3)=bbbr(12,i,j)
      bb(i3,j)=bb(i,j3)
      ab(i,j)=bbbr(8,i,j)
      ba(i,j)=ab(i,j)
      ab(i,j3)=bbbr(9,i,j)
      ba(i3,j)=ab(i,j3)
      ab(i3,j)=bbbr(10,i,j)
      ba(i,j3)=ab(i3,j)
      ab(i3,j3)=bbbr(11,i,j)
      ba(i3,j3)=ab(i3,j3)
160   continue
      call COUPL (aa,ab,ba,bb,cs,6,6,12,12)
      return
      end

#--------------------------------------------------------------------
      SUBROUTINE COULCF (bbb,qpi,a,xk,eps)
#
#     GENERAL CALCULATION OF
#     COULOMB-COEFFICIENTS BY THE EWALD METHOD
#     REF: MARADUDIN ET.AL., SOLID STATE PHYSICS, SUPP.3, 2ND.ED.,
#          "LATTICE DYNAMICS IN THE HARMONIC APPROXIMATION", CH.6.2
#          USING THE FORMULAS (6.2.18,21,22,24,53,64)
#
#     BBB ... COULOMB COEFFICIENT MATRIX
#             INCLUDING THE IRREGULAR PART|
#     QPI ... WAVEVECTOR MULTIPLIED BY PI.
#     A ..... TRANSLATIONAL BASIS VECTORS
#             1ST INDEX ARE COORDINATES, 2ND ARE ATOM #
#     XK  ... THE DISPLACEMENT VECTOR X(K')-X(K)
#     EPS ... CONVERGENCE PARAMETER FOR WIDTH OF GAUSSIAN
#             SEE MARADUDIN P.204
#
#--------------------------------------------------------------------
      real*4 erff
      real qpi(3), a(3,3), xk(3)
      complex bbb(3,3)
      dimension b(3,3),krm(3),kdm(3),tauk(3),xl(3),xlk(3),hab(3,3)
      complex phfac
      logical rzero
      data pi,pi2,pi4,pisqr/3.1415926535898,6.2831853072,12.5663706144,1
     1.7724538509/
#     CONVERGENCE PARAMETER FOR ACCURACY ABOUT 10**(-8)
      data convrg/4.17/
#
#     RECIPROCAL SPACE BASIS, VOLUME OF UNIT CELL
#
      do 100 i=1,3
      i1=mod(i,3)+1
      i2=mod(i+1,3)+1
      do 100 j=1,3
      j1=mod(j,3)+1
      j2=mod(j+1,3)+1
      b(i,j)=a(j1,i1)*a(j2,i2)-a(j1,i2)*a(j2,i1)
100   continue
      va=a(1,1)*b(1,1)+a(1,2)*b(2,1)+a(1,3)*b(3,1)
      do 110 i=1,3
      do 110 j=1,3
      b(i,j)=b(i,j)/va
110   continue
      va=abs(va)
#
#     CALCULATION OF PARAMETERS
#
#     SQRP**2 IS THE WIDTH OF THE GAUSSIAN (SEE MARADUDIN P.204)
      sqrp=eps/(va**.33333333)
      p4=4.0*sqrp**2
      pva=eps**3
      rzero=(xk(1)**2+xk(2)**2+xk(3)**2).lt.1.0e-10
#     AT QPI = 0 WE CHOOSE TO GO ALONG (100)-DIRECTION
      if (qpi(1)**2+qpi(2)**2+qpi(3)**2.lt.1.0e-10) qpi(1)=0.00001
      do 120 i=1,3
      do 120 j=1,3
      bbb(i,j)=(0.0,0.0)
120   continue
#
#     CALCULATION OF MAXIMAL LATTICE VECTORS
#
      do 130 i=1,3
      kdm(i)=ifix(convrg/sqrp*sqrt(b(i,1)**2+b(i,2)**2+b(i,3)**2)+0.5)
      krm(i)=ifix(convrg*sqrp/pi*sqrt(a(1,i)**2+a(2,i)**2+a(3,i)**2)+0.5
     1)
130   continue
#
#     SUMMATION IN RECIPROCAL SPACE
#     (MARADUDIN EQ.6.2.18)
#
      krm1=2*krm(1)+1
      krm2=2*krm(2)+1
      krm3=2*krm(3)+1
      do 160 k1=1,krm1
      fk1=float(k1-krm(1)-1)
      do 160 k2=1,krm2
      fk2=float(k2-krm(2)-1)
      do 160 k3=1,krm3
      fk3=float(k3-krm(3)-1)
      do 140 j=1,3
#     TAU + K
      tauk(j)=pi2*(fk1*b(1,j)+fk2*b(2,j)+fk3*b(3,j))+qpi(j)
140   continue
#     (TAU + K)*XK
      taukr=tauk(1)*xk(1)+tauk(2)*xk(2)+tauk(3)*xk(3)
#     (TAU + K)**2
      tauk2=tauk(1)**2+tauk(2)**2+tauk(3)**2
      tauk2e=tauk2/p4
      if (tauk2e.gt.convrg**2.or.tauk2e.eq.0.0) go to 160
      phfac=cexp(cmplx(-tauk2e,-taukr))
      do 150 i=1,3
      do 150 j=1,3
      bbb(i,j)=bbb(i,j)+cmplx(tauk(i)*tauk(j)/tauk2*pi4,0.0)*phfac
150   continue
160   continue
#
#     SUMMATION IN DIRECT SPACE
#     (MARADUDIN EQ.6.2.21)
#
      kdm1=2*kdm(1)+1
      kdm2=2*kdm(2)+1
      kdm3=2*kdm(3)+1
      do 200 k1=1,kdm1
      fk1=float(k1-kdm(1)-1)
      do 200 k2=1,kdm2
      fk2=float(k2-kdm(2)-1)
      do 200 k3=1,kdm3
      fk3=float(k3-kdm(3)-1)
      do 170 j=1,3
#     X(L')-X(L)
      xl(j)=fk1*a(j,1)+fk2*a(j,2)+fk3*a(j,3)
#     X(L'K')-X(LK)
      xlk(j)=xl(j)+xk(j)
170   continue
      z=xlk(1)**2+xlk(2)**2+xlk(3)**2
#     SQRT(P)*/X(L'K')-X(LK)/
      zr=sqrp*sqrt(z)
      if (zr.gt.convrg) go to 200
#     IF X(L'K') = X(LK), AND XK =0, SKIP THIS TERM
      if (rzero.and.abs(fk1)+abs(fk2)+abs(fk3).lt.1.0e-10) go to 200
#
#     CALCULATE H(ALFA,BETA)
#
      exq=exp(-zr**2)
      a1=(1.0-erff(zr)+2.0*zr/pisqr*exq)/zr**3*pva
      a2=3.0*a1+4.0/pisqr*exq*pva
      do 180 i=1,3
      do 180 j=1,3
      hab(i,j)=xlk(i)*xlk(j)/z*a2
#     EXTRA DIAGONAL TERM
      if (i.eq.j) hab(i,j)=hab(i,j)-a1
180   continue
#
      phfac=cexp(cmplx(.0,qpi(1)*xl(1)+qpi(2)*xl(2)+qpi(3)*xl(3)))
      do 190 i=1,3
      do 190 j=1,3
      bbb(i,j)=bbb(i,j)+phfac*cmplx(-hab(i,j),0.0)
190   continue
200   continue
#     FROM (6.2.53) WE GET AN EXTRA TERM FOR XK = 0
#     WHICH WAS EXCLUDED ABOVE
      if (.not.rzero) go to 220
      a1=-4.0/3.0/pisqr*pva
      do 210 i=1,3
      bbb(i,i)=bbb(i,i)+cmplx(a1,0.0)
210   continue
#
#     PHASE FACTOR
220   phfac=cexp(cmplx(0.0,qpi(1)*xk(1)+qpi(2)*xk(2)+qpi(3)*xk(3)))
      do 230 i=1,3
      do 230 j=1,3
      bbb(i,j)=bbb(i,j)*phfac
230   continue
      return
      end

#-----------------------------------------------------------------------
      SUBROUTINE COUPL (a11,a12,a21,a22,a,i1,j1,k1,l1)
#-----------------------------------------------------------------------
      real a11(i1,j1), a12(i1,j1), a21(i1,j1), a22(i1,j1), a(k1,l1)
#
#     PUT THE MATRICES A11,A12,A21,A22 INTO THE BIG MATRIX A
#                 (A11  A12)
#             A = (A21  A22)
#
      do 100 j=1,j1
      do 100 i=1,i1
      k=i+i1
      l=j+j1
      a(i,j)=a11(i,j)
      a(i,l)=a12(i,j)
      a(k,j)=a21(i,j)
      a(k,l)=a22(i,j)
100   continue
      return
      end

#-----------------------------------------------------------------------
      SUBROUTINE CSSCAL (n, sa, cx, incx)
#  Purpose:    Multiply a complex vector by a single-precision scalar,
#              y = ay.
#  Usage:      CALL CSSCAL (N, SA, CX, INCX)
#  Arguments:
#     N      - Length of vectors X.  (Input)
#     SA     - Real scalar.  (Input)
#     CX     - Complex vector of length MAX(N*IABS(INCX),1).
#                 (Input/Output)
#              CSSCAL replaces X(I) with SA*X(I) for I = 1,...,N.
#              X(I) refers to a specific element of CX.
#     INCX   - Displacement between elements of CX.  (Input)
#              X(I) is defined to be CX(1+(I-1)*INCX). INCX must be
#              greater than 0.
#-----------------------------------------------------------------------
#
      integer    n, incx
      real       sa
      complex    cx(*)
      integer    i, nincx
      intrinsic  cmplx
      complex    cmplx
#
      if (n .gt. 0) then
         if (incx .ne. 1) then
#                                  CODE FOR INCREMENT NOT EQUAL TO 1
            nincx = n*incx
            do 10  i=1, nincx, incx
               cx(i) = cmplx(sa,0.0e0)*cx(i)
   10       continue
         else
#                                  CODE FOR INCREMENT EQUAL TO 1
            do 20  i=1, n
               cx(i) = cmplx(sa,0.0e0)*cx(i)
   20       continue
         end if
      end if
      return
      end

#-----------------------------------------------------------------------
      SUBROUTINE MAT33 (a,b,c,d,e,f,g,h,p,q,s)
#     PUT ELEMENTS INTO A 3X3 MATRIX, AND MULTIPLY BY PHASE FACTOR
#-----------------------------------------------------------------------
#
#             (A  B  C)
#      S = Q* (D  E  F)
#             (G  H  P)
#
      real s(3,3)
      s(1,1)=q*a
      s(1,2)=q*b
      s(1,3)=q*c
      s(2,1)=q*d
      s(2,2)=q*e
      s(2,3)=q*f
      s(3,1)=q*g
      s(3,2)=q*h
      s(3,3)=q*p
      return
      end

#-----------------------------------------------------------------------
      SUBROUTINE MATOP (oper,a,na,ma,b,nb,mb,c,factor)
#
#     MATRIX MANIPULATIONS
#
#     PARAMETERS:
#     OPER ... VARIABLE OF TYPE "INTEGER" CONTROLLING THE ACTION
#              TAKEN BY MATOP:
#              '.    ' - B=A*FACTOR
#              '+    ' - C=A+B
#              '-    ' - C=A-B
#              '*    ' - C=A*B
#              'INV  ' - B=INVERSE OF A (USING IMSL)
#              'TRANS' - B=(TRANSPOSE OF A)*FACTOR
#              'EIG  ' - HERMITIAN EIGENVALUE PROBLEM (USING IMSL)
#              THE COMPLEX MATRIX (REAL,IMAG) IS (A,B)
#              THE EIGENVALUES ARE STORED IN THE DIAGONAL OF C
#              THE EIGENVECTORS (REAL,IMAG) ARE STORED IN (A,B)
#     A ...... MATRIX OF DIMENSION (NA,MA)
#     B ...... MATRIX OF DIMENSION (NB,MB)
#     C ...... MATRIX OF DIMENSION (NA,MB)
#     FACTOR . REAL NUMBER USED AS A MULTIPLICATIVE FACTOR
#              IN SOME OF THE OPERATIONS
#-----------------------------------------------------------------------
      character oper*5
      real a(na,ma),b(nb,mb),c(na,mb)

      real lambda(12)
#     COMPLEX V(6,6),AA(6,6)
      real v(12,12),aa(12,12)

# DETERMINE OPERATION TYPE
      if (oper.eq. '.    ') go to 110
      if (oper.eq. '+    '.or.oper.eq. '-    ') go to 150
      if (oper.eq. '*    ') go to 180
      if (oper.eq. 'INV  ') go to 210
      if (oper.eq. 'TRANS') go to 270
      if (oper.eq. 'EIG  ') go to 300
      write (6,100) oper  ! operation error
100   format (37h0improper operation in MATOP, oper = ,a5)
      stop

# *******  SCALAR MULTIPLICATION B = A*FACTOR  *************
# B MAY EQUAL A IN THIS CASE
# C IS A DUMMY ARRAY IN THIS CASE

110   if (na.eq.nb.and.ma.eq.mb) go to 130
      write (6,120) na,ma,nb,mb  ! dimension error
120   format (40h0dimension error in MATOP, na,ma,nb,mb =,4i10)
      stop
130   do 140 j=1,ma
      do 140 i=1,na
      b(i,j)=a(i,j)*factor
140   continue
      return

# ******* MATRIX ADDITION *******  
# C MAY EQUAL A OR B IN CALLING MATOP
# FACTOR IS A DUMMY

150   sign=1.0
      if (oper.eq. '-    ') sign=-1.0
      if (na.eq.nb.and.ma.eq.mb) go to 160
      write (6,120) na,ma,nb,mb   ! dimension error
      stop
160   do 170 j=1,ma
      do 170 i=1,na
      c(i,j)=a(i,j)+b(i,j)*sign
170   continue
      return
#
# ******* MATRIX MULTIPLICATION C = A*B  ******* 
# C MUST NOT EQUAL A OR B IN CALLING MATOP
# FACTOR IS A DUMMY

180   if (ma.eq.nb) go to 190
      write (6,120) na,ma,nb,mb  ! dimension error
      stop
190   do 200 j=1,mb
      do 200 i=1,na
      c(i,j)=0.0
      do 200 k=1,ma
      c(i,j)=c(i,j)+a(i,k)*b(k,j)
200   continue
      return

# ******* MATRIX INVERSION B = A**(-1)   USING LINRG (IMSL)  ******* 
# B MAY EQUAL A IN CALLING MATOP
# C AND FACTOR ARE DUMMIES

210   if (na.ne.ma.or.ma.ne.nb.or.nb.ne.mb) then
        write (6,120) na,ma,nb,mb     ! dimension error
        stop
      endif
220   if (na.gt.12) then
        write (6,230)
230     format (49h0na.gt.12 in MATOP, please increase dimension l,m)
        stop
      else
#        CALL LINRG(NA,A,NA,B,NA)
        call invmat(a,b,na,na)
      endif
      return

# ******* TRANSPOSE B = A(T)*FACTOR  ******* 
# B MAY NOT EQUAL A
# C IS A DUMMY

270   if (na.eq.mb.and.ma.eq.nb) go to 280
      write (6,120) na,ma,nb,mb  ! dimension error
      stop
280   do 290 j=1,ma
      do 290 i=1,na
      b(j,i)=a(i,j)*factor
290   continue
      return

# ******* HERMITIAN EIGENVALUE PROBLEM  ******* 
300   if (na.eq.ma.and.ma.eq.nb.and.nb.eq.mb) go to 310
      write (6,120) na,ma,nb,mb   ! dimension error
      stop
310   if (na.gt.6) then
        write (6,320)
320     format (43h0please increase dim lambda,vr,vi in MATOP/)
        stop
      else
        do j=1,na
          do i=1,na
            aa(i,j)=a(i,j)
            aa(i,j+na) = -b(i,j)
            aa(i+na,j) = b(i,j)
            aa(i+na,j+na) = a(i,j)
          end do
        end do
        na2 = 2*na
      
        call eigen(aa,na2,na2,lambda,v)

# ***** eigenvectors should already be normalized
#        DO J=1,6
#          WK1=1.0/SCNRM2(6,V(1,J),1)  ! normalize eigenvectors
#          CALL CSSCAL(6,WK1,V(1,J),1)
#        END DO
        do j=1,na
          c(j,j)=lambda(2*j-1)
          do i=1,na
            a(i,j)=v(i,2*j-1)
            b(i,j)=v(i+na,2*j-1)
          end do
        end do
      endif
      return
      end

#-----------------------------------------------------------------------
      SUBROUTINE PARBCM (itype,iprint)
#
#     INITIALIZE THE MODEL PARAMETERS
#
#     INPUT:
#     ITYPE ... 1-DIAMOND (default parameters)
#               2-SILICON  (default parameters)
#               3-GERMANIUM (natural, M=72.59, default parameters)
#               4-ALFA-TIN  (default parameters)
#               5-Germanium 70 (default parameters)
#               6-GE76  (default parameters)
#               11- Diamond with pars set in NEWPAR
#               12,13,14,15, etc
#     IPRINT .. 0-PARAMETERS NOT PRINTED
#               ELSE- PARAMETERS PRINTED ON TAPE6.
#
#     OUTPUT: 1/3 OF A PAGE ON TR12APE6, IF SPECIFIED BY IPRINT.
#-----------------------------------------------------------------------

      common /pram/ fcbc(6),fcion(6),mass,scale,acell2,z2eps,expram,thr
     1     am,celast(3,2)
      real mass


#     A ... PRIMITIVE TRANSLATIONAL VECTORS.
#           A(I,J) IS COORDINATE I OF THE J'TH VECTOR.
#     R12 . CONTAINS ALL THE RELATIVE COORDINATES OF THE IONS AND
#           BOND CHARGES, I.E. X(K')-X(K).
#           R12(I,J) IS COORDINATE J OF RELATIVE VECTOR I.

      common /basis/ a(3,3),r12(12,3)
      data a/1.,1.,0.,1.,0.,1.,0.,1.,1./  !1st 3 are 1st column
      data r12/ 0., .5, .25,-.25, .25,-.25, .5, .0, .5,-.5, .0, .5,
     1          0., .5, .25, .25,-.25,-.25, .0, .5, .5, .5, .5, .0,
     2          0., .5, .25,-.25,-.25, .25, .5, .5, .0, .0,-.5,-.5/

      data eesu/4.80325/,amu24/1.66053/
      data pi/3.1415926535898/

      itipe=mod(itype-1,10)+1      ! next 100 lines changed by STZO

      if (itipe.eq.1) then         ! DIAMOND PARAMETERS
        if (itype.eq.1) then
          fcbc(1)=51.682
          fcbc(2)=45.256
          fcbc(3)=4.669
          fcbc(4)=3.138
          fcbc(5)=-1.607
          fcbc(6)=3.138
          fcion(1)=-18.161
          fcion(2)=-6.02
          fcion(3)=0.0
          fcion(4)=0.0
          fcion(5)=0.0
          fcion(6)=0.0
          z2eps=0.885
        endif
        mass=12.0
        acell2=1.78
#     WARREN ET.AL., INELASTIC SCATTERING OF NEUTRONS, CONF. BOMBAY 1964
#                    IAEA VIENNA 1965, VOL.1,P.361
#    ALSO: BESERMAN (1972) (RAMAN, UNPUBLISHED)
        expram=39.96
#     MAC SKIMIN, PHYS.REV.105,116 (1957)
        celast(1,1)=10.76
        celast(2,1)=1.25
        celast(3,1)=5.76
      else if (itipe.eq.2) then   !  SILICON PARAMETERS
        if (itype.eq.2) then
          fcbc(1)=10.77
          fcbc(2)=2.17
          fcbc(3)=2.15
          fcbc(4)=2.15
          fcbc(5)=-2.15
          fcbc(6)=2.15
          fcion(1)=4.57
          fcion(2)=7.03
          fcion(3)=0.0
          fcion(4)=0.0
          fcion(5)=0.0
          fcion(6)=0.0
          z2eps=0.179
        endif
        mass=28.086
        acell2=2.715
#     DOLLING, INELASTIC SCATTERING OF NEUTRONS IN SOLIDS AND LIQUIDS,
#              IAEA VIENNA 1963, VOL.II, P.37
#              IAEA VIENNA 1965, VOL.I , P.249
        expram=15.53
#     MAC SKIMIN, PHYS.REV.105, 116 (1957)
#                 J.APPL.PHYS.24, 988 (1953)
        celast(1,1)=1.6563
        celast(2,1)=0.6390
        celast(3,1)=0.7956
      else if (itipe.eq.3) then   !  GERMANIUM PARAMETERS
        if (itype.eq.3) then
          fcbc(1)=9.9075
          fcbc(2)=1.5075
          fcbc(3)=2.10
          fcbc(4)=2.10
          fcbc(5)=-2.10
          fcbc(6)=2.10
          fcion(1)=5.1354
          fcion(2)=7.3512
          fcion(3)=0.0
          fcion(4)=0.0
          fcion(5)=0.0
          fcion(6)=0.0
          z2eps=0.162
        endif
        mass=72.59
        acell2=2.825
#     NILSSON, PHYS.REV.B3, 364 (1971)
        expram=9.11
#     MAC SKIMIN, J.APPL.PHYS.24, 988 (1953)
        celast(1,1)=1.2897
        celast(2,1)=0.48346
        celast(3,1)=0.67130
      else if (itipe.eq.4) then   !  ALFA-TIN PARAMETERS
        if (itype.eq.4) then
          fcbc(1)=9.49
          fcbc(2)=1.69
          fcbc(3)=1.95
          fcbc(4)=1.95
          fcbc(5)=-1.95
          fcbc(6)=1.95
          fcion(1)=5.94
          fcion(2)=8.18
          fcion(3)=0.0
          fcion(4)=0.0
          fcion(5)=0.0
          fcion(6)=0.0
          z2eps=0.163
        endif
        mass=118.69
        acell2=3.23
#     BUCHENAUER, PHYS.REV.B3, 1243 (1971)
        expram=5.97
#     PRICE, PHYS.REV.B1, 1268 (1971)
        celast(1,1)=0.690
        celast(2,1)=0.293
        celast(3,1)=0.362
      else if (itipe.eq.5) then   !  GERMANIUM 70 PARAMETERS
        if (itype.eq.5) then      !  added by STZO
          fcbc(1)=9.9075    ! These are just like natural Ge
          fcbc(2)=1.5075    !
          fcbc(3)=2.10      !
          fcbc(4)=2.10      !
          fcbc(5)=-2.10     !
          fcbc(6)=2.10      !
          fcion(1)=5.1354   !
          fcion(2)=7.3512   !
          fcion(3)=0.0      !
          fcion(4)=0.0      !
          fcion(5)=0.0      !
          fcion(6)=0.0      !
          z2eps=0.162       ! These are just like natural Ge
        endif
        mass=70.00    ! only the mass is different
        acell2=2.825
#     NILSSON, PHYS.REV.B3, 364 (1971)
        expram=9.11   ! natural Ge
#     MAC SKIMIN, J.APPL.PHYS.24, 988 (1953)
        celast(1,1)=1.2897       ! natural Ge
        celast(2,1)=0.48346      ! natural Ge
        celast(3,1)=0.67130      ! natural Ge
      else if (itipe.eq.6) then   !  GERMANIUM 76 PARAMETERS
        if (itype.eq.6) then      !  added by STZO
          fcbc(1)=9.9075    ! These are just like natural Ge
          fcbc(2)=1.5075    !
          fcbc(3)=2.10      !
          fcbc(4)=2.10      !
          fcbc(5)=-2.10     !
          fcbc(6)=2.10      !
          fcion(1)=5.1354   !
          fcion(2)=7.3512   !
          fcion(3)=0.0      !
          fcion(4)=0.0      !
          fcion(5)=0.0      !
          fcion(6)=0.0      !
          z2eps=0.162       ! These are just like natural Ge
        endif
        mass=75.69    ! only the mass is different
        acell2=2.825
#     NILSSON, PHYS.REV.B3, 364 (1971)
        expram=9.11   ! natural Ge
#     MAC SKIMIN, J.APPL.PHYS.24, 988 (1953)
        celast(1,1)=1.2897       ! natural Ge
        celast(2,1)=0.48346      ! natural Ge
        celast(3,1)=0.67130      ! natural Ge
      endif
#
#     CALCULATE DERIVED QUANTITIES
#
      fak=4.0*acell2**4/eesu**2
#     CONVERT INTO UNITS E**2/VCELL
      do 190 i=1,3
      celast(i,1)=celast(i,1)*fak
190   continue
#     SCALEFACTOR FOR FREQUENCY NU IN THZ (NOT OMEGA=2*PI*NU ||)
      scale=1.0e4/(2.0*pi)**2*eesu**2/((2.0*acell2)**3/4.0)/amu24
#     THEORETICAL ELASTIC CONSTANTS
#     Z2EPS = Z**2 / EPSILON
      p1=0.5*fcbc(1)+fcion(1)
      celast(1,2)=p1+2.*fcbc(3)+2.218*z2eps+8.*fcion(3)
      celast(2,2)=-p1+fcbc(2)+2.*fcion(2)+2.*fcbc(4)-fcbc(3)-fcbc(5)-30.
     1765*z2eps+8.0*fcion(4)-4.0*fcion(3)-4.0*fcion(5)
      celast(3,2)=p1-(fcbc(2)/2.+fcion(2)-17.653*z2eps)**2/p1+fcbc(3)+fc
     1bc(5)-1.109*z2eps+4.*fcion(3)+4.*fcion(5)
#     THEORETICAL RAMAN FREQUENCY
      thram=sqrt((4.0*fcbc(1)+8.0*fcion(1))*scale/mass)
      if (iprint.eq.0) return
#
#     PRINT OUT PARAMETERS
#
#     WRITE (6,200) ITIPE,FCBC,FCION,MASS,Z2EPS,ACELL2,
#    &              EXPRAM,THRAM,CELAST
200   format (1h , 'PARAMETERS FOR MATERIAL NO. ',i3,/
     &        1h , 'BOND CHARGE FORCES = ',6f12.4, ' E**2/VCELL',/,
     &        1h , 'ION FORCES = ',8x,6f12.4, ' E**2/VCELL',/,
     &        1h , 'MASS = ',f12.4, ' AMU',/,
     &        1h , 'Z**2/EPSILON= ',f12.4,/,
     &        1h , '1/2 CUBE LENGTH = ',f12.4, ' ANGSTROM',/
     &        1h , 'EXP.RAMAN FREQUENCY = ',f12.4, ' THZ',/,
     &        1h , 'TH. RAMAN FREQUENCY = ',f12.4, ' THZ,'/
     &        1h , 'EXP.ELASTIC CONSTANTS C11,C12,C44 = ',3f12.4,
     &                        'E**2/VCELL',/,
     &        1h , 'TH. ELASTIC CONSTANTS C11,C12,C44 = ',3f12.4,
     &                        'E**2/VCELL')
      return
      end

#-----------------------------------------------------------------------
      SUBROUTINE RMAT (rr,ri)
#
#     ION-ION FORCE DYNAMICAL MATRIX, DENOTED "R" IN W.WEBER'S ARTICLE
#     OUTPUT:
#     RR,RI ... REAL & IMAG PART OF MATRIX R
#-----------------------------------------------------------------------
      real rr(6,6), ri(6,6)

      common /pram/ fcbc(6),fcion(6),mass,scale,acell2,z2eps,expram,thr
     1am,celast(3,2)
      real mass
      common /qvect/ qpi(3)
      real mup4,nup4,lambp4
      cosx2=cos(qpi(1)/2.0)
      sinx2=sin(qpi(1)/2.0)
      cosy2=cos(qpi(2)/2.0)
      siny2=sin(qpi(2)/2.0)
      cosz2=cos(qpi(3)/2.0)
      sinz2=sin(qpi(3)/2.0)
      do 100 j=1,6
      do 100 i=1,6
      rr(i,j)=0.0
      ri(i,j)=0.0
100   continue
      alfa4=4.0*fcbc(1)
      alfap4=4.0*fcion(1)
      betap4=4.0*fcion(2)
      mup4=4.0*fcion(3)
      nup4=4.0*fcion(4)
      lambp4=4.0*fcion(5)
      deltp4=4.0*fcion(6)
#
#     CONSTRUCT MATRIX R FOR 1.N.N.
#
#     BY TRANSLATIONAL INVARIANCE WE HAVE:
      do 110 i=1,6
      rr(i,i)=alfap4+alfa4
110   continue
      rr(1,4)=-alfap4*cosx2*cosy2*cosz2
      rr(2,5)=rr(1,4)
      rr(3,6)=rr(1,4)
      rr(1,5)=betap4*sinx2*siny2*cosz2
      rr(1,6)=betap4*sinx2*cosy2*sinz2
      rr(2,6)=betap4*cosx2*siny2*sinz2
      rr(2,4)=rr(1,5)
      rr(3,4)=rr(1,6)
      rr(3,5)=rr(2,6)
      ri(1,4)=alfap4*sinx2*siny2*sinz2
      ri(2,5)=ri(1,4)
      ri(3,6)=ri(1,4)
      ri(1,5)=-betap4*cosx2*cosy2*sinz2
      ri(1,6)=-betap4*cosx2*siny2*cosz2
      ri(2,6)=-betap4*sinx2*cosy2*cosz2
      ri(2,4)=ri(1,5)
      ri(3,4)=ri(1,6)
      ri(3,5)=ri(2,6)
#
#     CONSTRUCT MATRIX R FOR 2.N.N.
#
      cosx=cos(qpi(1))
      sinx=sin(qpi(1))
      cosy=cos(qpi(2))
      siny=sin(qpi(2))
      cosz=cos(qpi(3))
      sinz=sin(qpi(3))
      rr2=mup4*(2.0-cosx*(cosy+cosz))+lambp4*(1.0-cosy*cosz)
      rr(1,1)=rr(1,1)+rr2
      rr(4,4)=rr(4,4)+rr2
      rr2=mup4*(2.0-cosy*(cosz+cosx))+lambp4*(1.0-cosz*cosx)
      rr(2,2)=rr(2,2)+rr2
      rr(5,5)=rr(5,5)+rr2
      rr2=mup4*(2.0-cosz*(cosx+cosy))+lambp4*(1.0-cosx*cosy)
      rr(3,3)=rr(3,3)+rr2
      rr(6,6)=rr(6,6)+rr2
      rr2=nup4*sinx*siny
      rr(1,2)=rr(1,2)+rr2
      rr(4,5)=rr(4,5)+rr2
      rr2=nup4*sinx*sinz
      rr(1,3)=rr(1,3)+rr2
      rr(4,6)=rr(4,6)+rr2
      rr2=nup4*siny*sinz
      rr(2,3)=rr(2,3)+rr2
      rr(5,6)=rr(5,6)+rr2
      ri(1,2)=-deltp4*sinz*(cosx-cosy)
      ri(1,3)=-deltp4*siny*(cosx-cosz)
      ri(2,3)=-deltp4*sinx*(cosy-cosz)
      ri(4,5)=-ri(1,2)
      ri(4,6)=-ri(1,3)
      ri(5,6)=-ri(2,3)
      do 120 i=1,5
      i1=i+1
      do 120 j=i1,6
      ri(j,i)=-ri(i,j)
      rr(j,i)=rr(i,j)
120   continue
      return
      end

#-----------------------------------------------------------------------
      SUBROUTINE SMATX (s)
#
#     BC-BC FORCE DYNAMICAL MATRIX, DENOTED "S" IN W.WEBER'S ARTICLE
#     OUTPUT:
#     S ... THE MATRIX S (A REAL MATRIX)
#-----------------------------------------------------------------------
      real s(12,12)

      common /pram/ fcbc(6),fcion(6),mass,scale,acell2,z2eps,expram,thr
     1am,celast(3,2)
      real mass
      common /qvect/ qpi(3)
      real a11(3,3), a12(3,3), a21(3,3), a22(3,3), b11(6,6), b12(6,
     16), b21(6,6), b22(6,6)
      a=2.0*fcbc(1)+4.0*fcbc(3)+2.0*fcbc(5)
      b=2.0*fcbc(2)+2.0*fcbc(4)
      c=fcbc(3)
      d=fcbc(4)
      e=fcbc(5)
      f=fcbc(6)
      call MAT33 (a,b,b,b,a,b,b,b,a,1.0,a11)
      call MAT33 (a,-b,b,-b,a,-b,b,-b,a,1.0,a22)
      q=-2.0*cos((qpi(1)+qpi(3))/2.0)
      call MAT33 (c,f,d,-f,e,-f,d,f,c,q,a12)
      do 100 i=1,3
      do 100 j=1,3
      a21(j,i)=a12(i,j)
100   continue
      call COUPL (a11,a12,a21,a22,b11,3,3,6,6)
      call MAT33 (a,-b,-b,-b,a,b,-b,b,a,1.,a11)
      call MAT33 (a,b,-b,b,a,-b,-b,-b,a,1.,a22)
      q=-2.0*cos((qpi(1)-qpi(3))/2.0)
      call MAT33 (c,-f,-d,f,e,-f,-d,f,c,q,a12)
      do 110 i=1,3
      do 110 j=1,3
      a21(j,i)=a12(i,j)
110   continue
      call COUPL (a11,a12,a21,a22,b22,3,3,6,6)
      q=-2.0*cos((qpi(2)+qpi(3))/2.0)
      call MAT33 (e,-f,-f,f,c,d,f,d,c,q,a11)
      q=-2.0*cos((qpi(1)+qpi(2))/2.0)
      call MAT33 (c,d,f,d,c,f,-f,-f,e,q,a12)
      q=-2.0*cos((qpi(1)-qpi(2))/2.0)
      call MAT33 (c,-d,f,-d,c,-f,-f,f,e,q,a21)
      q=-2.0*cos((qpi(2)-qpi(3))/2.0)
      call MAT33 (e,f,-f,-f,c,-d,f,-d,c,q,a22)
      call COUPL (a11,a12,a21,a22,b12,3,3,6,6)
      do 120 i=1,6
      do 120 j=1,6
      b21(j,i)=b12(i,j)
120   continue
      call COUPL (b11,b12,b21,b22,s,6,6,12,12)
      return
      end

#-----------------------------------------------------------------------
      SUBROUTINE TMAT (tr,ti)
      real tr(6,12),ti(6,12)
#
#     ION-BC FORCE DYNAMICAL MATRIX, DENOTED "T" IN W.WEBER'S ARTICLE
#     OUTPUT:
#     TR,TI ... REAL & IMAG PART OF MATRIX T
#-----------------------------------------------------------------------

      common /pram/ fcbc(6),fcion(6),mass,scale,acell2,z2eps,expram,thr
     1am,celast(3,2)
      real mass
      common /qvect/ qpi(3)
      real a11(3,3), a12(3,3), a21(3,3), a22(3,3), tr1(6,6), tr2(6,
     16)
      a=fcbc(1)
      b=fcbc(2)
      q=-cos((qpi(1)+qpi(2)+qpi(3))/4.0)
      call MAT33 (a,b,b,b,a,b,b,b,a,q,a11)
      q=-cos((qpi(1)-qpi(2)+qpi(3))/4.0)
      call MAT33 (a,-b,b,-b,a,-b,b,-b,a,q,a12)
      call COUPL (a11,a12,a11,a12,tr1,3,3,6,6)
      q=-cos((qpi(1)-qpi(2)-qpi(3))/4.0)
      call MAT33 (a,-b,-b,-b,a,b,-b,b,a,q,a11)
      q=-cos((qpi(1)+qpi(2)-qpi(3))/4.0)
      call MAT33 (a,b,-b,b,a,-b,-b,-b,a,q,a12)
      call COUPL (a11,a12,a11,a12,tr2,3,3,6,6)
      do 100 j=1,6
      do 100 i=1,6
      tr(i,j)=tr1(i,j)
      tr(i,j+6)=tr2(i,j)
100   continue
      q=-sin((qpi(1)+qpi(2)+qpi(3))/4.0)
      call MAT33 (a,b,b,b,a,b,b,b,a,q,a11)
      call MAT33 (a,b,b,b,a,b,b,b,a,-q,a21)
      q=sin((qpi(1)-qpi(2)+qpi(3))/4.0)
      call MAT33 (a,-b,b,-b,a,-b,b,-b,a,q,a12)
      call MAT33 (a,-b,b,-b,a,-b,b,-b,a,-q,a22)
      call COUPL (a11,a12,a21,a22,tr1,3,3,6,6)
      q=-sin((qpi(1)-qpi(2)-qpi(3))/4.0)
      call MAT33 (a,-b,-b,-b,a,b,-b,b,a,q,a11)
      call MAT33 (a,-b,-b,-b,a,b,-b,b,a,-q,a21)
      q=sin((qpi(1)+qpi(2)-qpi(3))/4.0)
      call MAT33 (a,b,-b,b,a,-b,-b,-b,a,q,a12)
      call MAT33 (a,b,-b,b,a,-b,-b,-b,a,-q,a22)
      call COUPL (a11,a12,a21,a22,tr2,3,3,6,6)
      do 110 j=1,6
      do 110 i=1,6
      ti(i,j)=tr1(i,j)
      ti(i,j+6)=tr2(i,j)
110   continue
      return
      end


#-----------------------------------------------------------------------
      SUBROUTINE eigen(a,n,np,eival,eivect)
#-----------------------------------------------------------------------

      real a(np,np),eival(np),eivect(np,np)
      real d(12),e(12)
#
      do j=1,np
        do i=1,np
          eivect(i,j) = a(i,j)
        enddo
      enddo
#
      call tred2(eivect,n,np,d,e)
      call tqli(d,e,n,np,eivect)
      call eigsrt(d,eivect,n,np)
#
      do i=1,np
        eival(i) = d(i)
      enddo
#
      return
      end


#-----------------------------------------------------------------------
      SUBROUTINE TRED2(a,n,np,d,e)
#-----------------------------------------------------------------------

      real a(np,np),d(np),e(np)
#
      if(n.gt.1)then
        do 18 i=n,2,-1
          l=i-1
          h=0.
          scale=0.
          if(l.gt.1)then
            do 11 k=1,l
              scale=scale+abs(a(i,k))
11          continue
            if(scale.eq.0.)then
              e(i)=a(i,l)
            else
              do 12 k=1,l
                a(i,k)=a(i,k)/scale
                h=h+a(i,k)**2
12            continue
              f=a(i,l)
              g=-sign(sqrt(h),f)
              e(i)=scale*g
              h=h-f*g
              a(i,l)=f-g
              f=0.
              do 15 j=1,l
                a(j,i)=a(i,j)/h
                g=0.
                do 13 k=1,j
                  g=g+a(j,k)*a(i,k)
13              continue
                if(l.gt.j)then
                  do 14 k=j+1,l
                    g=g+a(k,j)*a(i,k)
14                continue
                endif
                e(j)=g/h
                f=f+e(j)*a(i,j)
15            continue
              hh=f/(h+h)
              do 17 j=1,l
                f=a(i,j)
                g=e(j)-hh*f
                e(j)=g
                do 16 k=1,j
                  a(j,k)=a(j,k)-f*e(k)-g*a(i,k)
16              continue
17            continue
            endif
          else
            e(i)=a(i,l)
          endif
          d(i)=h
18      continue
      endif
      d(1)=0.
      e(1)=0.
      do 23 i=1,n
        l=i-1
        if(d(i).ne.0.)then
          do 21 j=1,l
            g=0.
            do 19 k=1,l
              g=g+a(i,k)*a(k,j)
19          continue
            do 20 k=1,l
              a(k,j)=a(k,j)-g*a(k,i)
20          continue
21        continue
        endif
        d(i)=a(i,i)
        a(i,i)=1.
        if(l.ge.1)then
          do 22 j=1,l
            a(i,j)=0.
            a(j,i)=0.
22        continue
        endif
23    continue
      return
      end


#-----------------------------------------------------------------------
      SUBROUTINE TQLI(d,e,n,np,z)
#-----------------------------------------------------------------------

      real d(np),e(np),z(np,np)
#
      if (n.gt.1) then

        do 11 i=2,n

          e(i-1)=e(i)
11      continue
        e(n)=0.
        do 15 l=1,n
          iter=0
1         do 12 m=l,n-1
            dd=abs(d(m))+abs(d(m+1))
            if (abs(e(m))+dd.eq.dd) go to 2
12        continue
          m=n
2         if(m.ne.l)then
           
            if(mod(iter+1,100).eq.0) then                    
#
#
             write(*,*)   'TQLI warning: too many iterations '                    
          endif   
            iter=iter+1
            g=(d(l+1)-d(l))/(2.*e(l))
            r=sqrt(g**2+1.)
            g=d(m)-d(l)+e(l)/(g+sign(r,g))
            s=1.
            c=1.
            p=0.
            do 14 i=m-1,l,-1
              f=s*e(i)
              b=c*e(i)
              if(abs(f).ge.abs(g))then
                c=g/f
                r=sqrt(c**2+1.)
                e(i+1)=f*r
                s=1./r
                c=c*s
              else
                s=f/g
                r=sqrt(s**2+1.)
                e(i+1)=g*r
                c=1./r
                s=s*c
              endif
              g=d(i+1)-p
              r=(d(i)-g)*s+2.*c*b
              p=s*r
              d(i+1)=g+p
              g=c*r-b
              do 13 k=1,n
                f=z(k,i+1)
                z(k,i+1)=s*z(k,i)+c*f
                z(k,i)=c*z(k,i)-s*f
13            continue
14          continue
            d(l)=d(l)-p
            e(l)=g
            e(m)=0.
            go to 1
          endif
15      continue
      endif
      return
      end

#-----------------------------------------------------------------------
      SUBROUTINE EIGSRT(d,v,n,np)
#-----------------------------------------------------------------------
#
      dimension d(np),v(np,np)
#
      do 13 i=1,n-1
        k=i
        p=d(i)
        do 11 j=i+1,n
          if(d(j).ge.p)then
            k=j
            p=d(j)
          endif
11      continue
        if(k.ne.i)then
          d(k)=d(i)
          d(i)=p
          do 12 j=1,n
            p=v(j,i)
            v(j,i)=v(j,k)
            v(j,k)=p
12        continue
        endif
13    continue
      return
      end


#-----------------------------------------------------------------------
      real*4 FUNCTION MYERF(x)
# renamed by J.S. => avoid conflict with some intrinsic libraries 
#-----------------------------------------------------------------------
      if(x.lt.0.)then
        MYERF=-GAMMP(.5,x**2)
      else
        MYERF=GAMMP(.5,x**2)
      endif
      return
      end
#
#-----------------------------------------------------------------------
      FUNCTION GAMMP(a,x)
#-----------------------------------------------------------------------
      if(x.lt.0..or.a.le.0.)pause
      if(x.lt.a+1.)then
        call GSER(gamser,a,x,gln)
        GAMMP=gamser
      else
        call GCF(gammcf,a,x,gln)
        GAMMP=1.-gammcf
      endif
      return
      end

#-----------------------------------------------------------------------
      SUBROUTINE GCF(gammcf,a,x,gln)
#-----------------------------------------------------------------------
      parameter (itmax=100,eps=3.e-7)
#
      gln=GAMMLN(a)
      gold=0.
      a0=1.
      a1=x
      b0=0.
      b1=1.
      fac=1.
      do 11 n=1,itmax
        an=float(n)
        ana=an-a
        a0=(a1+a0*ana)*fac
        b0=(b1+b0*ana)*fac
        anf=an*fac
        a1=x*a0+anf*a1
        b1=x*b0+anf*b1
        if(a1.ne.0.)then
          fac=1./a1
          g=b1*fac
          if(abs((g-gold)/g).lt.eps)go to 1
          gold=g
        endif
11    continue
      pause  'A too large, ITMAX too small'
1     gammcf=exp(-x+a*alog(x)-gln)*g
      return
      end

#-----------------------------------------------------------------------
      FUNCTION GAMMLN(xx)
#-----------------------------------------------------------------------
      real*8 cof(6),stp,half,one,fpf,x,tmp,ser
      data cof,stp/76.18009173d0,-86.50532033d0,24.01409822d0,
     *    -1.231739516d0,.120858003d-2,-.536382d-5,2.50662827465d0/
      data half,one,fpf/0.5d0,1.0d0,5.5d0/
#
      x=xx-one
      tmp=x+fpf
      tmp=(x+half)*log(tmp)-tmp
      ser=one
      do 11 j=1,6
        x=x+one
        ser=ser+cof(j)/x
11    continue
      GAMMLN=tmp+log(stp*ser)
      return
      end

#-----------------------------------------------------------------------
      SUBROUTINE GSER(gamser,a,x,gln)
#-----------------------------------------------------------------------
      parameter (itmax=100,eps=3.e-7)
#
      gln=GAMMLN(a)
      if(x.le.0.)then
        if(x.lt.0.)pause
        gamser=0.
        return
      endif
      ap=a
      SUM=1./a
      del=SUM
      do 11 n=1,itmax
        ap=ap+1.
        del=del*x/ap
        SUM=SUM+del
        if(abs(del).lt.abs(SUM)*eps)go to 1
11    continue
      pause  'A too large, ITMAX too small'
1     gamser=SUM*exp(-x+a*log(x)-gln)
      return
      end


#-----------------------------------------------------------------------
      SUBROUTINE invmat(a,b,n,np)
#-----------------------------------------------------------------------
      real a(np,np),b(np,np)
      real c(20,20)
      integer indx(12)
#
      do j=1,n
        do i=1,n
          b(i,j) = .0
          c(i,j) = a(i,j)
        enddo
        b(j,j) = 1.
      enddo
#
      call ludcmp(a,n,np,indx,d)
#
      do j=1,np
        call lubksb(a,n,np,indx,b(1,j))
      enddo
#
      do j=1,n
        do i=1,n
          a(i,j) = c(i,j)
        enddo
      enddo
      return
      end

#-----------------------------------------------------------------------
       SUBROUTINE LUDCMP(a,n,np,indx,d)
#-----------------------------------------------------------------------
      parameter (nmax=100,tiny=1.0e-20)
      dimension a(np,np),indx(n),vv(nmax)
#
      d=1.
      do 12 i=1,n
        aamax=0.
        do 11 j=1,n
          if (abs(a(i,j)).gt.aamax) aamax=abs(a(i,j))
11      continue
        if (aamax.eq.0.) pause  'Singular matrix.'
        vv(i)=1./aamax
12    continue
      do 19 j=1,n
        if (j.gt.1) then
          do 14 i=1,j-1
            SUM=a(i,j)
            if (i.gt.1)then
              do 13 k=1,i-1
                SUM=SUM-a(i,k)*a(k,j)
13            continue
              a(i,j)=SUM
            endif
14        continue
        endif
        aamax=0.
        do 16 i=j,n
          SUM=a(i,j)
          if (j.gt.1)then
            do 15 k=1,j-1
              SUM=SUM-a(i,k)*a(k,j)
15          continue
            a(i,j)=SUM
          endif
          dum=vv(i)*abs(SUM)
          if (dum.ge.aamax) then
            imax=i
            aamax=dum
          endif
16      continue
        if (j.ne.imax)then
          do 17 k=1,n
            dum=a(imax,k)
            a(imax,k)=a(j,k)
            a(j,k)=dum
17        continue
          d=-d
          vv(imax)=vv(j)
        endif
        indx(j)=imax
        if(j.ne.n)then
          if(a(j,j).eq.0.)a(j,j)=tiny
          dum=1./a(j,j)
          do 18 i=j+1,n
            a(i,j)=a(i,j)*dum
18        continue
        endif
19    continue
      if(a(n,n).eq.0.)a(n,n)=tiny
      return
      end

#-----------------------------------------------------------------------
      SUBROUTINE LUBKSB(a,n,np,indx,b)
#-----------------------------------------------------------------------
      dimension a(np,np),indx(n),b(n)
#
      ii=0
      do 12 i=1,n
        ll=indx(i)
        SUM=b(ll)
        b(ll)=b(i)
        if (ii.ne.0)then
          do 11 j=ii,i-1
            SUM=SUM-a(i,j)*b(j)
11        continue
        else if (SUM.ne.0.) then
          ii=i
        endif
        b(i)=SUM
12    continue
      do 14 i=n,1,-1
        SUM=b(i)
        if(i.lt.n)then
          do 13 j=i+1,n
            SUM=SUM-a(i,j)*b(j)
13        continue
        endif
        b(i)=SUM/a(i,i)
14    continue
      return
      end


#-----------------------------------------------------------------------
      real*4 FUNCTION erff(x)
# fast version of erf(x), using lookup tables
# added by J.S.
#-----------------------------------------------------------------------
      implicit none
      integer*4 mt
      real*4 xmax,dx,x
      parameter (mt=1024,xmax=5.0)
      parameter (dx=xmax/mt)
      real*4 erftab(0:mt),a,b,xx,ss,z
      real*4 myerf
      integer*4 init,nx,i
      data init/0/
      save erftab
      
      if (init.eq.0) then
        erftab(0)=0.0
      do i=1,mt
          erftab(i)=myerf(i*dx)
      enddo
      init=1
      endif
      
      xx=abs(x)
      if (xx.le.xmax) then
        ss=sign(1.0,x)
      z=xx/dx
        nx=nint(z)
      if (nx.lt.1) nx=nx+1
      if (nx.ge.mt) nx=nx-1         
        a=(erftab(nx+1)+erftab(nx-1))/2.0-erftab(nx)
      b=(erftab(nx+1)-erftab(nx+1))/2.0
      erff=ss*(erftab(nx)+a*(z-nx)**2+b*(z-nx))          
      else
        erff=1.0      
      endif

      return
      end
#
# $Log: res_exci_bcm.f,v $
# Revision 1.10  2006/05/11 16:33:10  saroun
# debugging bcm and phon models
#
# Revision 1.9  2006/05/10 18:46:32  saroun
# *** empty log message ***
#
# Revision 1.8  2006/05/08 21:39:57  saroun
# x   change BCM model: optional width profile as damped oscillators
# x   BUG FIX: CFG prompt does not show current as default
# x   BUG FIX: po navratu z MFIT, MAPSQ se vypne ray-tracing (MFIT->FIT)
# x   BUG FIX: subsequent print shows ellipsoid instead of  cloud
# x   PLOT jenom resol. funkce + disp. surface (bez dat)
# x   BUG FIX: SQOm, kdyz neni zadan komentar, zobrazi se divne znaky
#