src/sim_cal.f

Fortran project SIMRES, source module src/sim_cal.f.

Source module last modified on Mon, 28 Mar 2005, 19:05;
HTML image of Fortran source automatically generated by for2html on Mon, 23 May 2005, 21:29.


#//////////////////////////////////////////////////////////////////////
#////
#////  R E S T R A X   4.1
#////
#////  Subroutines from the RESCAL program
#////  + some transformation procedures
#////
#////  * SUBROUTINE RECLAT
#////  * FUNCTION   TRANS(A,I)
#////  * SUBROUTINE CONVRT(A,AZ)
#////  * SUBROUTINE TRANSMAT
#////  * SUBROUTINE QNORM(Q,QRLU,QANG)
#////  * SUBROUTINE TRNVCT
#////  
#//////////////////////////////////////////////////////////////////////

#****************************
        SUBROUTINE RECLAT
#****************************
#        11-05-78
#        25-06-79 UPDATED
#        OVERLAY #11
#        COMPUTATION OF THE TRANSFORMATION MATRIX B
#        AND ORIENTATION MATRIX S
#        REQUIRED: VECTOR A=UNIT CELL SIDES, VECTOR ALFA=CELL ANGLES
#        A1=VECTOR IN DIRECTION ACO, A2=VECTOR IN DIRECTION BCO
#        X=SUM(S(1,I)*Q(I)) AND Y=SUM(S(2,I)*Q(I)) ARE CARTESIAN
#        COORDINATES IN SCATTERING PLANE, TRANSFORMED FROM VECTOR
#        Q W.R.T. CRYSTAL LATTICE.  SUM(S(3,I)*Q(I)) SHOULD=0
#
#       S(i,j) transforms Q(hkl) in the r.l.u. into 
#       the orthogonal coordinates given by A(i),B(i) vectors 
#       defining the scattering plane
#
#       SD(i,j) is a matrix defined so that
#       Qhkl = SQRT((hkl)*SD*(hkl))   (in A^-1)
#
#***********************************************************************
      implicit none     
      INCLUDE 'const.inc'
      INCLUDE 'inout.inc'
      INCLUDE 'rescal.inc'
      real*8 eps   
      parameter(eps=1.0d-8)
      real*8 s(3,3),sd(3,3),b(3),cosb(3)
      common /s/s,sd,b,cosb
      integer*4 ier
      common /error/ier
      real*8 aq(3),alfa(3),a1(3),a2(3)

      real*8 v1(3),v2(3),v3(3),u(3,3),cosa(3),sina(3),sinb(3)
      real*8 a(3),bb(3,3),zd,rd,cc,c1,c2,c3,ss
      integer*4 i,j,k,l,m

      do i=1,3
        aq(i)=res_dat(i_as+i-1)
        alfa(i)=res_dat(i_aa+i-1)
        a1(i)=res_dat(i_ax+i-1)
        a2(i)=res_dat(i_bx+i-1)
      enddo
        zd=2.d0*pi
        rd=2*pi/360.d0
        cc=0.
        do 1 i=1,3
          a(i)=aq(i)/zd
          if(abs(a(i)).gt..00001) goto 6
          
#       write(*,*)I,A(I)   
          ier=1
          go to 100
6          cosa(i)=cos(alfa(i)*rd)
          sina(i)=sin(alfa(i)*rd)
    1          cc=cc+cosa(i)*cosa(i)
        cc=1.+2.*cosa(1)*cosa(2)*cosa(3)-cc
        if(cc.le.0.) goto 23
        cc=sqrt(cc)
        j=2
        k=3
        do 2 i=1,3
        b(i)=sina(i)/(a(i)*cc)             ! length of the r.l. axes
        cosb(i)=(cosa(j)*cosa(k)-cosa(i))/(sina(j)*sina(k))
        sinb(i)=sqrt(1.-cosb(i)*cosb(i))   ! angles btw. r.l. axes
        j=k
    2        k=i
    
    
        do i=1,3   
           sd(i,i)=b(i)**2
        end do
        sd(1,2)=(cosa(1)*cosa(2)-cosa(3))/aq(1)/aq(2)*(2*pi/cc)**2
        sd(1,3)=(cosa(1)*cosa(3)-cosa(2))/aq(1)/aq(3)*(2*pi/cc)**2
        sd(2,3)=(cosa(2)*cosa(3)-cosa(1))/aq(2)/aq(3)*(2*pi/cc)**2
        sd(2,1)=sd(1,2)
        sd(3,1)=sd(1,3)
        sd(3,2)=sd(2,3)    
    
    
        bb(1,1)=b(1)
        bb(2,1)=0.
        bb(3,1)=0.
        bb(1,2)=b(2)*cosb(3)
        bb(2,2)=b(2)*sinb(3)
        bb(3,2)=0.
        bb(1,3)=b(3)*cosb(2)
        bb(2,3)=-b(3)*sinb(2)*cosa(1)
        bb(3,3)=1/a(3)                     ! checked ... O.K.
#        GENERATION OF ORIENTATION MATRIX REC. LATTICE TO SCATTERING PLANE
        do 3 i=1,3
        c1=0.
        c2=0.
        do 4 j=1,3
        c1=c1+bb(i,j)*a1(j)
    4        c2=c2+bb(i,j)*a2(j)
        v1(i)=c1
    3        v2(i)=c2
        v3(1)=v1(2)*v2(3)-v1(3)*v2(2)
        v3(2)=v1(3)*v2(1)-v1(1)*v2(3)
        v3(3)=v1(1)*v2(2)-v1(2)*v2(1)
        v2(1)=v3(2)*v1(3)-v3(3)*v1(2)
        v2(2)=v3(3)*v1(1)-v3(1)*v1(3)
        v2(3)=v3(1)*v1(2)-v3(2)*v1(1)
        c1=(v1(1)*v1(1)+v1(2)*v1(2)+v1(3)*v1(3))
        c2=(v2(1)*v2(1)+v2(2)*v2(2)+v2(3)*v2(3))
        c3=v3(1)**2+v3(2)**2+v3(3)**2
        if(abs(c1)-eps)14,14,15
   14        ier=5
        goto 100
   15        if(abs(c2)-eps)16,16,17
   16        ier=6
        goto 100
   17        if(abs(c3)-eps)18,18,19
   18        ier=7
        goto 100
   19        continue
        c1=sqrt(c1)
        c2=sqrt(c2)
        c3=sqrt(c3)
        do 5 i=1,3
        u(1,i)=v1(i)/c1
        u(2,i)=v2(i)/c2
    5        u(3,i)=v3(i)/c3
        do 7 k=1,3
        do 7 m=1,3
        ss=0.
        do 8 l=1,3
    8        ss=ss+u(k,l)*bb(l,m)
        s(k,m)=ss
    7        continue
  100        if(ier)200,200,20
   20        goto(22,22,22,23,24,24,24),ier
   22        write(sout,30)
        goto 200
   23        write(sout,31)
        ier=1
        goto 200
   24        write(sout,32)
  200        continue
        return
   30        format( ' RECLAT: Check Lattice Spacings (AS,BS,CS)'/)
   31        format( ' RECLAT: Check Cell Angles (AA,BB,CC)'/)
   32        format( ' RECLAT: Check Scattering Plane (AX....BZ)'/)
        end

        real*8 FUNCTION TRANS(a,i)
#***********************************************************************
#***********************************************************************
        implicit none
        integer*4 i
        real*8 a(3)
        real*8 s(3,3),sd(3,3),b(3),cosb(3)
        common /s/s,sd,b,cosb
        TRANS=s(i,1)*a(1)+s(i,2)*a(2)+s(i,3)*a(3)
        return
        end

#--------------------------------------------------------------------
      SUBROUTINE TRANSMAT
#  creates transformation matrices for the four coordinate systems:     
#  R  ...  reciprocal lattice  
#  C  ...  Cooper & Nathans coordinates (C&N)
#  G  ...  C&N rotated so that X(1) // grad E(Q)
#  D  ...  G rotated so that X(4) // normal to the disp. surface    
#--------------------------------------------------------------------
      implicit none
      INCLUDE 'const.inc'
      INCLUDE 'inout.inc'
      INCLUDE 'rescal.inc'
      real*8 s(3,3),sd(3,3),b(3),cosb(3)
      common /s/s,sd,b,cosb
      real*8 q3(3),z,xg,yg,zg,xq,yq,gg,gg1,qq,co,si,co1,si1
      real*8 gnr,gna,zra,phi,si2,co2,mgd(4,4)
      real*8 TRANS
      integer*4 i,j
      
#         write(*,*) 'TRANSMAT: ',Q
      z=sqrt(qhkl(1)**2+qhkl(2)**2+qhkl(3)**2)
      if(z.eq.0) then 
        q3(1)=1
        q3(2)=0
        q3(3)=0
      else
        do i=1,3
          q3(i)=qhkl(i)
        enddo
      endif            
      xg=TRANS(grd,1)
      yg=TRANS(grd,2)
      zg=TRANS(grd,3)
      xq=TRANS(q3,1)
      yq=TRANS(q3,2)                              
      gg=sqrt(xg**2+yg**2+zg**2)
      gg1=sqrt(xg**2+yg**2)
      qq=sqrt(xq**2+yq**2)
      
      co=xq/qq
      si=yq/qq
      do 11 i=1,3
        mcr(1,i)=s(1,i)*co+s(2,i)*si
        mcr(2,i)=s(2,i)*co-s(1,i)*si
        mcr(3,i)=s(3,i)
        mcr(4,i)=0.
11        mcr(i,4)=0.
      mcr(4,4)=1.                            ! (MCR): r.l.u. --> C&N
                  

#///  Transformation from (grad(E)//x) to C&N
      si1=zg/gg
      co1=sqrt(1-si1**2)
      si2=(xq*yg-yq*xg)/qq/gg1
      co2=(xg*xq+yg*yq)/qq/gg1
 
      mcg(1,1)= co1*co2
      mcg(2,1)= co1*si2
      mcg(3,1)= +si1     
      mcg(1,2)= -si2
      mcg(2,2)= co2
      mcg(3,2)=0 
      mcg(1,3)=-si1*co2 
      mcg(2,3)=-si1*si2 
      mcg(3,3)=co1
      do 10 i=1,3
        mcg(i,4)=0
        mcg(4,i)=0
10    continue
      mcg(4,4)=1                         ! (MCG): grad   -->  C&N
      
      do 20 i=1,4
        do 21 j=1,4
21        mgd(i,j)=0
        mgd(i,i)=1
20    continue

      call QNORM(grd(1),gnr,gna)
      zra=gnr/gna      
      
      phi=atan(res_dat(i_gmod)*zra)               ! units of GMOD are Energy/r.l.u.!

      mgd(1,1)=cos(phi)
      mgd(4,4)=cos(phi)
      mgd(1,4)=-sin(phi)
      mgd(4,1)=sin(phi)                  ! (MGD): grad   <--  disp 
         
      call MXM(1,4,4,mcg,mgd,mcd)        ! (MCD): C&N    <--  disp     
      call INVERT(4,mcr,4,mrc,4)         ! (MRC): C&N    -->  r.l.u.
      call MXM(-1,4,4,mcd,mcr,mdr)       ! (MDR): r.l.u. -->  disp

      return
      end


#--------------------------------------------
      SUBROUTINE QNORM(q,qrlu,qang)
#     returns norm of Qhkl in r.l.u and A^-1
#--------------------------------------------
      implicit none
      INCLUDE 'const.inc'
      INCLUDE 'inout.inc'
      INCLUDE 'rescal.inc'
      real*8 q(3),qrlu,qang
      real*8 s(3,3),sd(3,3),b(3),cosb(3)
      common /s/s,sd,b,cosb
      real*8 v3(3),V3xV3,TRANS
      integer*4 i    
      
      qrlu=sqrt(V3XV3(q,q)+2*q(1)*q(2)*cosb(3)+2*q(2)*q(3)*cosb(1)+
     1    2*q(1)*q(3)*cosb(2))                   ! norm of G(3) in r.l.u.
      
      do i=1,3
        v3(i)=TRANS(q,i)
      end do  
      qang= sqrt(v3(1)**2+v3(2)**2+v3(3)**2)    ! norm of G(3) in [A-1]
      
      return
      end  
      
#
#***********************************************************************
      SUBROUTINE RESMAT
#    calculates resolution matrices their parameters from MC 
#***********************************************************************
      implicit none
      INCLUDE 'const.inc'
      INCLUDE 'inout.inc'
      INCLUDE 'rescal.inc'
      real*8 DETERM
      
      integer*4 i,j,k,ncnt,icnt,ialloc
      real*8 pp, e(4),aux(4,4),suma,z1,z2
      real*8 cov(4,4),covrl(4,4),mean(4)
      real*8 volresmc,resmatcn(4,4),resmatrl(4,4)
      real*8 bragcn(4),vancn(4),bragrl(4),vanrl(4)
      common /resmatpar/ volresmc,resmatcn,resmatrl,bragcn,vancn,
     *                   bragrl,vanrl

1     format(a,i,6(2x,g11.5))

      do j=1,4
        do k=1,4
          resmatcn(j,k)=0.d0
          resmatrl(j,k)=0.d0
          cov(j,k)=0.d0
        enddo
        bragcn(j)=0.d0  
        bragrl(j)=0.d0  
        vancn(j)=0.d0  
        vanrl(j)=0.d0 
        mean(j)=0.d0
      enddo 
      ncnt=0
#// get number of accumulated (Q,E) events
      call NSTORE_N(i,ncnt,ialloc)
#      write(*,*) NCNT,IALLOC
#      pause
      if(ncnt.lt.3.or.ialloc.lt.ncnt) goto 99
      icnt=0
      suma=0.d0
#// accumulate covariance matrix
      do  i=1,ncnt
         call NSTORE_GETQE(i,e,pp,0.)
         do j=1,4
           do k=1,4
            cov(j,k)=cov(j,k)+e(j)*e(k)*pp
           enddo
           mean(j)=mean(j)+e(j)*pp
         enddo  
         suma=suma+pp
         icnt=icnt+1
      enddo 
      if (suma.le.0.or.icnt.le.1) goto 99  
#// normalize covariance matrix
      do j=1,4
         mean(j)=mean(j)/suma
         do k=1,4
           cov(j,k)=cov(j,k)/suma
         enddo
      enddo   
#// subtract mean 
      do j=1,4
         do k=1,4
           cov(j,k)=cov(j,k)-mean(j)*mean(k)
         enddo
      enddo   
#// invert COV to get resolution matrix in C&N coord.
      call INVERT(4,cov,4,resmatcn,4)
#// transform to reciprocal lattice coord.      
      call MXM(1,4,4,resmatcn,mcr,aux)
      call MXM(-1,4,4,mcr,aux,resmatrl) 
#// Bragg widths are calculated from diagonal elements of resol. mat.
      do j=1,4
        z1=abs(resmatcn(j,j))
        z2=abs(resmatrl(j,j))
        if (z1.gt.1.d-16) bragcn(j)=r8ln2/sqrt(z1)
        if (z2.gt.1.d-16) bragrl(j)=r8ln2/sqrt(z2)
      enddo
#// Vanad widths are calculated from diagonal elements of covar. mat.
      call INVERT(4,resmatrl,4,covrl,4)
      do j=1,4
        vancn(j)=r8ln2*sqrt(abs(cov(j,j)))
        vanrl(j)=r8ln2*sqrt(abs(covrl(j,j)))
      enddo
      volresmc=sqrt(abs(DETERM(cov,4,aux)))*4*pi**2
      return
      
99    continue
      
      end