exci/res_exci_phon.f

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

Source module last modified on Fri, 12 May 2006, 12:08;
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:
#-----------------------------------
# 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)
# Damped oscillators (up to 6 branches) with free dispersion gradient
# Phonon polarization is included in the calculation of S(Q,omega)
# Energy dispersion E(Q) as a quadratic form (in 3 dimensions)
# For wen(j)=0, zero intrinsic width is considered.
#------------------------------------------------------------------------------
      implicit none
      
#----------------------- *** DON'T CHANGE *** ------------------------------
      INCLUDE 'const.inc'
      INCLUDE 'exci.inc'
      integer*4 icom,excinit   
      real*8 q(4),omexc(6),sqom(6)
      
#-------------------------- *** END *** ------------------------------------
     
# **** Local user declarations ****

      integer*4 j,i
      real*8 qphon(3),qphon0(3),sig1(3),sig2(3),sig3(3),dqphi(3)
      real*8 dq1,dq2,dq3,bf,kt,z,f(6)
      character*1 ch(6) 

#// internal variables describing the model
#// This common is not shared with the rest of RESTRAX, only within this file. 
#// These variables are actually shared with the READEXCIPAR subroutine below.
      real*8 temp,tau(3),qnom(3),evec(3,6)
      real*8 fqsq(6),omega0(6),domdq(6),b1(6),b2(6),b3(6)
      common /excipar/ temp,tau,qnom,fqsq,evec,omega0,domdq,b1,b2,b3

# function defined in reclat.f, makes dot product of vectors in rec. lat. coordinates:
# REAL*8 QxQ(A,B) for REAL*8 A(3),B(3)
      real*8 QxQ      

#----------------------- *** DATA section *** ----------------------------     
      
# **** DEFAULT values of internal model variables ****
      data temp/1.5/
      data fqsq/6*1./
      data omega0/7.,15.,4*1.0/, domdq/0.,10.,4*0/
      data b1/-10,-10.,4*0/, b2/6*0/, b3/6*0/ 
      data ch/ '1', '2', '3', '4', '5', '6'/
      data tau/0.,0.,2/
      data qnom/0.1,0.1,2.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

        call SETEXCIDEFAULT   ! DONīT CHANGE

# Set model identification string: 
      phontitle= 'Phonons, quadratic form dispersion, damped oscillators'

#// Define fixed parameters (=0), default=1 (all free)
        fixparam(1)=0   ! let Intensity fixed !!

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

#// Initial widths in energy (default=1meV), 
#// Set wen(i)=0 for zero-width branches
#        WEN(1)=0.D0 
#        WEN(2)=2.D0

#**** 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 (default=exc.par)
        phonname= 'phon.par'
      
# Define names of free parameters for i>2:
#        PARNAME(3)='Position'
#        PARNAME(4)='Width'       

#// Number of free model parameters
#        NTERM=4
                
        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  
#
#// actual number of free model parameters
      nterm=2+nbr*4
#// Assign free model parameters to param():

      do i=1,nbr
        param(3+4*(i-1)) = fqsq(i)
        param(4+4*(i-1)) = omega0(i)
        param(5+4*(i-1)) = wen(i)
        param(6+4*(i-1)) = domdq(i)
      enddo
      
#// set parameter names for all branches
      do i=1,nbr
        parname(3+4*(i-1))= 'Int'//ch(i)
        parname(4+4*(i-1))= 'Omega'//ch(i)
        parname(5+4*(i-1))= 'Gamma'//ch(i)
        parname(6+4*(i-1))= 'Grad'//ch(i)
      enddo
      do i=nbr+1,6
        parname(3+4*(i-1))= ' '
        parname(4+4*(i-1))= ' '
        parname(5+4*(i-1))= ' '
        parname(6+4*(i-1))= ' '
      enddo
            
#----------------------- *** DON'T CHANGE *** -------------------------

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

#********************************************************************************
#                                                                   
#                   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.
# 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
      do i=1,nbr
        fqsq(i)=param(3+4*(i-1))
        omega0(i)=param(4+4*(i-1))
        wen(i)=param(5+4*(i-1))
        domdq(i)=param(6+4*(i-1))
      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 !! -------------------------

      do i=1,3
         qphon(i) = q(i)-tau(i)      ! actual phonon q
       qphon0(i) = qnom(i)-tau(i)  ! nominal phonon q
         dqphi(i) = q(i)-qnom(i)     ! difference from nominal phonon q
      enddo

# Use POLVECT to get polarization vectors:
# SIG1 ...  along qphon0
# SIG2 ...  perpendicular to qphon0, in plane
# SIG3 ...  perpendicular to qphon0, off plane (vertical)
      call POLVECT(qnom,tau,sig1,sig2,sig3,0) 

#* get projections of dqphi along qphon0, perpendicular to qphon0, and vertical
      dq1=QxQ(dqphi,sig1)
      dq2=QxQ(dqphi,sig2)
      dq3=QxQ(dqphi,sig3)
      
# include phonon polarization in S(Q,E)
      do j=1,nbr
        f(j) = 0.
        do i=1,3
            f(j)=f(j)+q(i)*evec(i,j)
        enddo
        f(j)=fqsq(j)*f(j)**2
      enddo

#// Define dispersion as a quadratic form:
#// domdq ... gradient along q
#// b1,b2,b3 ... curvartures along q, perpendicular and vertical, resp.
      do j=1,nbr
      omexc(j) = omega0(j)+dq1*domdq(j)+
     *           dq1**2*b1(j)+dq2**2*b2(j)+dq3**2*b3(j)
      enddo

#// Define scattering cross-section: damped oscillator
      kt=temp/11.609  ! 11.609 ... conversion kT -> meV
      z=exp(-q(4)/kt)  
      if (abs(1.0-z).lt.1.d-5) then
          bf=kt
      else
          bf=q(4)/(1-z)
      endif          
      do j=1,nbr
#// bf = Bose factor x omega
        if(fqsq(j).eq.0) then
            sqom(j)=0.
        else
          if(wen(j).ne.0) then           
              sqom(j) = f(j)*abs(wen(j))*bf/
     &        ((q(4)**2-omexc(j)**2)**2+wen(j)**2*q(4)**2)
          else
              if(abs(omexc(j)).lt.0.0001) omexc(j)=0.0001
              sqom(j) = bf*f(j)/omexc(j)**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'
      
      integer*4 i
      real*8 temp,tau(3),qnom(3),evec(3,6)
      real*8 fqsq(6),omega0(6),domdq(6),b1(6),b2(6),b3(6)
      common /excipar/ temp,tau,qnom,fqsq,evec,omega0,domdq,b1,b2,b3
      
      write(*,*)  'EXCI: '//phontitle

#// Report some model values:  
14    format( 'Temperature [K]: ',g10.4)
      write(*,14) temp
10    format( 'Number of branches: ',i2)
      write(*,10) nbr
11    format( 'tau (B.Z. center): ',3(1x,i3))
      write(*,11) (nint(tau(i)),i=1,3)
12    format( 'q0 (center of quadratic disp. surface): ',3(1x,g10.4))
      write(*,12) (qnom(i),i=1,3)
13    format( 'disp. curvatures (L, T, vertical): ',3(1x,g10.4))
      do i=1,nbr
        write(*,13) b1(i),b2(i),b3(i)
      enddo
      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'
          
      integer*4 i,j
      real*8 z,dum
      real*8 temp,tau(3),qnom(3),evec(3,6)
      real*8 fqsq(6),omega0(6),domdq(6),b1(6),b2(6),b3(6)
      common /excipar/ temp,tau,qnom,fqsq,evec,omega0,domdq,b1,b2,b3

      rewind(excunit)  ! call rewind for compatibility with g77
# read model parameters from file 
      read (excunit,*,err=998) nbr  ! number of branches (1..6)
      if (nbr.gt.6) nbr=6
      if (nbr.le.0) nbr=1
      read (excunit,*,err=998) temp   ! temperature
      read (excunit,*,err=998) tau    ! tau
      read (excunit,*,err=998) qnom   ! Q0, reference point for disp. branches
      do i=1,3
        tau(i)=1.d0*nint(tau(i))
      enddo 
      do i=1,nbr                         ! F(Q)^2, Polarization, Gamma, Omega0, Gradient, Curvatures
       read (excunit,*,err=998) fqsq(i),(evec(j,i),j=1,3),wen(i),
     &       omega0(i),domdq(i),b1(i),b2(i),b3(i)
      enddo      
      
# normalize polarization vectors:
      do j=1,nbr
        call QNORM(evec(1,j),z,dum)      
        do i=1,3
         evec(i,j)=evec(i,j)/z
        enddo
      enddo
     
      return

998   write(*,*)  'Format error?! Cannot read excitation parameters.'       
      return
      end