src/ness/ness_mos.f

Fortran project SIMRES, source module src/ness/ness_mos.f.

Source module last modified on Tue, 12 Apr 2005, 16:04;
HTML image of Fortran source automatically generated by for2html on Mon, 23 May 2005, 21:29.


#      PROGRAM MOSAIC
#      CALL TEST      
#      END

      SUBROUTINE TEST
      
      implicit none
      
      INCLUDE 'const.inc'
      INCLUDE 'structures.inc'  
      record /CRYSTAL/ cr
      real*8 GETMI,sqrt2pi
      parameter (sqrt2pi=2.506628275)

      real*8 lambda,eta,thick,chi,r0,fwhm,integ,weff,wmos,GETQKIN,reff
      real*8 ki,qkin,mi,temp
      integer*4 ierr
      character*8 crname,sdum
      data crname/ '        '/
      save crname



1     format( 'Crystal name [',a8, ']: ',$) 
2     format(a8)
3     format( 'ETA[' ']  D[mm]  chi[deg]  T[K]: ',$)          
4     format( 'ki = ',$)
5     format( 2x, 'ki          ',
     *         'fwhm        ',
     *         'R           ',
     *         'Integral    ',
     *         'Qkin.       ',
     *         'mi [cm^-1]  ',
     *         'sigma [b]   ')
           
           
15    format(10(f10.6,2x))
6     format( 'ETA[' ']  D[mm]  chi[deg]  T[K]: ',4(2x,f10.6))     
      
50    write(*,*)
      write(*,1) crname
      read(*,2) sdum
      call MakeUpCase(sdum)
      if ((sdum(1:3).eq. 'END').or.(sdum(1:3).eq. 'QUI')) then
         goto 200
      else if (sdum.ne. '        ') then
         crname=sdum                     
         call READCRYST(cr,crname)
      endif
      write(*,3)
      read(*,*) eta,thick,chi,temp
      cr.hmos=eta*pi/180/60/r8ln2
      cr.frame.size(3)=thick
      cr.chi=chi*pi/180
      ierr=0
      write(*,*) 
      write(*,*) crname
      write(*,6) eta,thick,chi,temp
      write(*,5) 
      do while (ierr.eq.0)
         write(*,4)
         read(*,*,err=100,iostat=ierr) ki
         if (ki.eq.0) goto 200
         lambda=2*pi/ki
       qkin= GETQKIN(cr,lambda)
       mi= GETMI(cr,lambda,temp)   
         call getrefpar(cr,lambda,qkin,mi,r0,integ,wmos,fwhm)
       if (r0.gt.0) then 
            weff=integ/sqrt2pi/r0*r8ln2
            reff=integ/sqrt2pi/wmos*r8ln2
         else
            weff=0
            reff=0
         endif
         weff=weff*180*60/pi
         wmos=wmos*180*60/pi
       fwhm=fwhm*180*60/pi
         write(*,15) ki,fwhm,r0,reff,qkin,mi,mi*cr.vol
      enddo
100   continue        
      goto 50
200   continue      
      end

#xxxxxxxxxCxxxxxxxxxCxxxxxxxxxCxxxxxxxxxCxxxxxxxxxCxxxxxxxxxCxxxxxxxxxCX
#--------------------------------------------------------
      real*8 FUNCTION GETQKIN(cr,lambda)
#--------------------------------------------------------
      implicit none
      INCLUDE 'const.inc'
      INCLUDE 'structures.inc'  
      record /CRYSTAL/ cr
      real*8 sthb,lambda,cthb
      
      sthb=lambda/2/cr.dhkl
      if (sthb.ge.9.99999d-1) then
         GETQKIN=0
       return
      endif    
      cthb=sqrt(1.d0-sthb**2)
      GETQKIN=cr.qml*cr.dhkl*sthb**2/pi/cthb
            
      end

#--------------------------------------------------------
      real*8 FUNCTION GETREFDYN(cr,lambda)
#--------------------------------------------------------
      implicit none
      INCLUDE 'const.inc'
      INCLUDE 'structures.inc'  
      record /CRYSTAL/ cr
      real*8 z,z1,z2,grad
      integer*4 i,j
      real*8 sthb,lambda,cthb,thb,qk,sigma,ki(3),kg(3)
      real*8 GETDW
      
      sthb=lambda/2/cr.dhkl      
      if (sthb.ge.9.99999d-1) then
         GETREFDYN=0
       return
      endif
                
      cthb=sqrt(1.d0-sthb**2)
      qk=cr.qml*cr.dhkl*sthb**2/pi   !  = Qkin*cos(thb)
      qk=qk*GETDW(cr,300.d0,z1,z2,lambda)/10.
      thb=asin(sthb)
      ki(1)=-cos(thb-cr.chi)
      ki(2)=0.
      ki(3)=-sin(thb-cr.chi)
      kg(1)=-cos(thb+cr.chi)
      kg(2)=0.
      kg(3)=+sin(thb+cr.chi)
      grad=0.
      do i=1,3
        z=0.
      do j=1,3
        z=z+cr.dg_dr(i,j)*ki(j)
      enddo
      grad=grad+z*kg(i)
      enddo
      grad=abs(grad)/cr.gtot 
      if (grad/qk.lt.1d-2) then
         z=1.d0
      else    
         sigma=qk/grad
         z=1.d0-exp(-sigma)
      endif
      GETREFDYN=z
      end

#-------------------------------------------------------------
      real*8 FUNCTION GETDTH(cr)
# Get the spread of Bragg angles due to the elastic deformation      
#-------------------------------------------------------------
      implicit none
      INCLUDE 'structures.inc'  
      record /CRYSTAL/ cr
      real*8 z,grad
      integer*4 i,j
      real*8 ki(3),kg(3),r(3),tin,tout,sthb,cthb    
      
      sthb=sin(abs(cr.thb))      
      if (sthb.ge.9.99999d-1) then
         GETDTH=0.d0
       return
      endif
          
      ki(1)=-cos(cr.thb-cr.chi)
      ki(2)=0.
      ki(3)=-sin(cr.thb-cr.chi)
      kg(1)=-cos(cr.thb+cr.chi)
      kg(2)=0.
      kg(3)=+sin(cr.thb+cr.chi)
      grad=0.
      do i=1,3
        z=0.
      do j=1,3
        z=z+cr.dg_dr(i,j)*ki(j)
      enddo
      grad=grad+z*kg(i)
      enddo
#      write(*,*) CR.FRAME.NAME,grad,cos(CR.THB),CR.GTOT
      cthb=sqrt(1.d0-sthb**2)
      if (cr.gtot.ge.1.d-10) then
         grad=abs(grad)/cthb/cr.gtot  
      else
         grad=0.d0
      endif
      r(1)=0.
      r(2)=0.
      r(3)=0.
      call CR_BORDER(cr,r,ki,tin,tout)
#      write(*,*) TIN,TOUT
      GETDTH=grad*(tout-tin)
      end
      
#--------------------------------------------------------------------
      real*8 FUNCTION GETEFFMOS(cr)
# Get effective mosaicity 
# = change of refl. plane orientation along incident beam path
#--------------------------------------------------------------------
      implicit none
      INCLUDE 'structures.inc'  
      record /CRYSTAL/ cr
      real*8 z,g01,g02
      integer*4 i,j
      real*8 ri(3),rf(3),g1(3),g2(3),cthb,schi,cchi  
      
      cthb=cos(cr.thb-cr.chi)    
      ri(1)=0.5*cthb*cr.frame.size(3)
      ri(2)=0.
      ri(3)=0.5*cr.frame.size(3)
      rf(1)=-0.5*cthb*cr.frame.size(3)
      rf(2)=0.
      rf(3)=-0.5*cr.frame.size(3)
      do i=1,3
        g1(i)=0.
        g2(i)=0.
      do j=1,3
        g1(i)=g1(i)+cr.dg_dr(i,j)*ri(j)
        g2(i)=g2(i)+cr.dg_dr(i,j)*rf(j)
      enddo
      enddo
      schi=sin(cr.chi)
      cchi=sqrt(1.d0-schi**2)
      g1(1)=g1(1)+cr.gtot*schi
      g1(3)=g1(3)-cr.gtot*cchi
      g2(1)=g2(1)+cr.gtot*schi
      g2(3)=g2(3)-cr.gtot*cchi
      z=0.
      do i=1,3
        z=z+g1(i)*g2(i)
      enddo
      if (abs(g01*g02).lt.1.d-10) then
         z=0.d0
      else
        g01=sqrt(g1(1)**2+g1(3)**2)
        g02=sqrt(g2(1)**2+g2(3)**2)
        z=acos(z/g01/g02)
      endif
      GETEFFMOS=z
      end
      
#--------------------------------------------------------
      real*8 FUNCTION GETDW(cr,t,b0,bt,lambda)
#///  Debye-Waller factor
#/// cubic crystals, Debye model      
#--------------------------------------------------------
      implicit none
      INCLUDE 'const.inc'
      INCLUDE 'structures.inc'  
      record /CRYSTAL/ cr
      integer*4 p
      real*8 b0,bt,x,z,dx,ksi,t,e,lambda
      integer*4 i
      b0=2873/cr.a/cr.thetad
      x=cr.thetad/t
      z=0
      dx=x/100
      do i=1,101
        if (i.eq.1.or.i.eq.101) then
          p=1
        else if (mod(i,2).eq.1) then
          p=4
        else
          p=2
        endif
        if (i.eq.1) then
          z=-1
        else  
          ksi=(i-1)*dx
          z=z+p*ksi/(exp(ksi)-1)
        endif
      enddo
      z=z*dx/3
      bt=4*b0*z/x**2
      e=hsqov2m*(2*pi/lambda)**2
#      GETDW=exp(-(B0+BT)*CR.C2*E)
      GETDW=exp(-2.0*(b0+bt)/(cr.dhkl*2)**2)  
      
      end
      
#--------------------------------------------------------
      real*8 FUNCTION GETMI(cr,lambda,t)
#//  calculates absorption according to A.Freund 
#//  (Nucl.Inst.Meth. 213 (1983) 495-501)
#
#//  sig1 .... multiple-phonon
#//  sig2 .... single-phonon      

#//  sigmaa ... absorption for 1A neutrons [barn*A^-1]
#//  sigmab ... bound-atom scattering cross-section [barn]
#//  V0 .... volume [A^3]/atom
#//  A  .... atomic number
#//  thetaD .... Debye temperature (K)
#//  C2 .... constant from the Freund's paper  [A^-2 meV^-1]
#--------------------------------------------------------
      implicit none
      INCLUDE 'const.inc'
      INCLUDE 'structures.inc'  
      record /CRYSTAL/ cr
      real*8 kb
      parameter (kb  = 0.086174)    ! Boltzman constant [meV*K^-1]     
      real*8 lambda,e,t,r,x,ifact,b(0:30)
      real*8 sig0,sig1,sig2,GETDW,b0,bt,sigtot
      integer*4 i

      data b /1,-0.5,0.166667,0,-0.033333,0,0.0238095,0,-0.033333,
     *  0,0.0757576,0,-0.253114,0,1.16667,0,-7.09216,0,54.9712,
     *  0,-529.124,0,6192.12,0,-86580.3,0,1.42551717e6,0,
     * -2.7298231e7,0,6.01580874e8/
      
      GETMI=0
      e=hsqov2m*(2*pi/lambda)**2
      
      
      x=cr.thetad/t
      r=0
      ifact=1
      do i=0,20
        ifact=ifact*i
      if (i.eq.0) ifact=1
        r=r+b(i)*x**(i-1)/(ifact*(i+5/2))
      enddo
      r=3*(r+b(22)*x**21/(ifact*21*22*(22+5/2))/2)          
      
      
      cr.dw=GETDW(cr,t,b0,bt,lambda)
      sig0=cr.sigmaa*lambda          
      sig1=cr.sigmab*(cr.a/(cr.a+1))**2*(1-exp(-(b0+bt)*cr.c2*e))
      if (x.le.6) then
        sig2=cr.sigmab/cr.a*sqrt(kb*cr.thetad/e)*r
      else
        sig2=cr.sigmab/cr.a*sqrt(kb*cr.thetad/e)*3.3/sqrt(x**7)
      endif
      sigtot=(sig0+sig1+sig2)/cr.vol
      GETMI=sigtot
#      IF(INDEX(CR.FRAME.NAME,'PG').GT.0) THEN  ! Pyrolytic graphite
#         IF(sigtot*CR.VOL.lt.4) GETMI=4/CR.VOL
#      ENDIF
      
#10    FORMAT(a13,3(G13.5,2x))
#      write(*,10) 'sigma     : ',sig0,sig1,sig2   
#      write(*,10) 'DW,B0,BT  : ',DW,B0,BT   
#      write(*,10) 'R,t/E,s/A : ',R,KB*CR.thetaD/E,CR.sigmab/CR.A   
#      write(*,10) 'V0, MI    : ',CR.VOL,(sig0+sig1+sig2)/CR.VOL 
      
      end
     
#-----------------------------------------------------------------
      SUBROUTINE GETREFPAR(cr,lambda,qkin,mi0,r0,integ,wg,fwhm)
#-----------------------------------------------------------------
      implicit none
      INCLUDE 'const.inc'
      INCLUDE 'structures.inc'  
      record /CRYSTAL/ cr
      real*8 lambda,r0,fwhm,integ,x(101),y(101),sumsq,MOSREF
      real*8 mi0,qkin,maxy,x2,x1,dx,wg,xmin,sthb,thb,cthb,stmchi,stpchi
      integer*4 i,i1,i2
      
      if(qkin.le.0) then
        r0=0
      integ=0
      wg=0
      fwhm=0
      return
      endif         

      maxy=0
      if(cr.hmos.gt.sec) then
        dx=cr.hmos*10./101.
        xmin=-cr.hmos*5
      else   
         sthb=lambda/2/cr.dhkl
         if (abs(sthb).ge.1) then
           return
         endif
         thb=asin(sthb)
         cthb=cos(thb)
         stmchi=abs(sin(thb-cr.chi))
         stpchi=abs(sin(thb+cr.chi))
         dx=cr.rh*cos(cr.chi)/cthb*(1-(1+cr.poi)*stpchi*stmchi)
       dx=dx+cr.dgr*1.d-4*cos(cr.dga-cr.chi)/cthb
       if(stmchi.gt.0.05) then
          dx=dx*cr.frame.size(3)/stmchi
       else
          dx=dx*cr.frame.size(1)/stpchi
       endif   
       dx=abs(dx)*1.5
       xmin=-dx/2
       dx=dx/101
      endif    

      do i=1,101
        x(i)=xmin+dx*i
        y(i)=MOSREF(cr,x(i),lambda,mi0,qkin)
      if (maxy.lt.y(i)) maxy=y(i)
      enddo
      i=1
      do while ((i.le.101).and.(y(i).le.maxy/2))
         i=i+1
      enddo
      i1=i
      i=101
      do while ((i.ge.1).and.(y(i).le.maxy/2))
         i=i-1
      enddo
      i2=i      
      x1=x(i1-1)+(maxy/2-y(i1-1))/(y(i1)-y(i1-1))*dx
      x2=x(i2)+(maxy/2-y(i2))/(y(i2+1)-y(i2))*dx
      fwhm=x2-x1
          
       
      
      integ=y(1)+y(101)
      do i=2,100
         if (mod(i,2).eq.0) then
            integ=integ+4*y(i)
         else
            integ=integ+2*y(i)
         endif
      enddo
      integ=integ/3*(x(2)-x(1))
      sumsq=y(1)*x(1)**2+y(101)*x(101)**2
      do i=2,100
         if (mod(i,2).eq.0) then
            sumsq=sumsq+4*y(i)*x(i)**2
         else
            sumsq=sumsq+2*y(i)*x(i)**2
         endif
      enddo
      sumsq=sumsq/3*(x(2)-x(1))
      wg=sqrt(sumsq/integ)*r8ln2
      r0=MOSREF(cr,0.d0,lambda,mi0,qkin)
      end
                    


#-------------------------------------------------
      SUBROUTINE READCRYST(cr,crname)
#// QML ... Maier-Leibnitz reflectivity in A^-1 cm^-1
#// QML = 4*PI*(F*dhkl/V0)**2       
#-------------------------------------------------
      implicit none
      INCLUDE 'inout.inc'  
      INCLUDE 'structures.inc'  
      record /CRYSTAL/ cr
      character*8 crname1
      character*(*) crname
      character*128 line
      integer*4 indx,ires,is,il,iline
1     format(a)      
2     format( 'Format error in the crystal library, line ',i5 )      
3     format( 'Crystal ',a, ' not found in the library' )      
4     format( 'Cannot open crystal library' )      
      indx=0
      
      call BOUNDS(crname,is,il)
      crname1=crname(is:is+il-1)
      call MakeUpCase(crname1)        
            
      call OPENRESFILE( 'crystal.lib',1,ires,1) 
      if(ires.le.0) goto 99
#      OPEN(unit=1,file='crystal.lib',status='OLD',err=98)
      iline=0
      do while (indx.eq.0)
         read(1,1,end=100,err=97) line
         iline=iline+1
         call MakeUpCase(line)
         indx=index(line,crname1)
      enddo
100   close(1)      
      if(indx.eq.0) goto 98
      read(line(indx+8:len(line)),*,err=97) 
     * cr.dhkl,cr.qml,cr.sigmab,cr.sigmaa,cr.vol,cr.a,cr.thetad,cr.c2,
     * cr.poi
      cr.c2=cr.c2/1000   ! from eV^-1 to meV^-1
      return

# format error      
97    write(smes,2) iline      
      cr.vol=0
      close(1)      
      return

# crystal not found
98    write(smes,3) crname1(1:il)
      cr.vol=0 
      close(1)      
      return
            
# cannot open      
99    write(smes,4)
      cr.vol=0    
      return
      

      end 
      

#-------------------------------------------------
      real*8 FUNCTION MOSREF(cr,x,lambda,mi0,qkin)
#-------------------------------------------------
      implicit none
      INCLUDE 'const.inc'
      INCLUDE 'structures.inc'  
      real*8 sqrt2pi
      parameter(sqrt2pi=2.506628275)
      real*8 x,lambda
      record /CRYSTAL/ cr
      real*8 eta1,stpchi,stmchi,thb,sthb,cthb,mi0,cgeom
      real*8 gamma1,gamma2,qkin,sigma,rb,rl,dum
      
      
      eta1=cr.hmos  
      sthb=lambda/2/cr.dhkl
      if (abs(sthb).ge.1) then
         MOSREF=0
         return
      endif
      thb=asin(sthb)
      cthb=cos(thb)
      stmchi=abs(sin(thb-cr.chi))
      stpchi=abs(sin(thb+cr.chi))
      
      if (eta1.gt.sec.and.(stmchi.lt.1.d-5.or.
     *    stpchi.lt.1.d-5)) then
        write(*,*)  'Cannot calculate reflectivity '
        MOSREF=0
        return
      endif        
      if (eta1.le.sec) then      ! bent perfect crystal 
         dum=cr.rh*cos(cr.chi)/cthb*(1-(1+cr.poi)*stpchi*stmchi)
       dum=dum+cr.dgr*1.d-4*cos(cr.dga-cr.chi)/cthb    
       if(stmchi.gt.0.05) then
          dum=dum*cr.frame.size(3)/stmchi
       else
          dum=dum*cr.frame.size(1)/stpchi
       endif   
         if (abs(dum)/2.gt.abs(x)) then
          MOSREF=cr.ref
       else
          MOSREF=0
       endif
      else            ! mosaic crystal
#!!! sign of CHI is oposite to the definition in TRAX
        cgeom=(stpchi-stmchi)/(stmchi+stpchi)
        gamma1=(1/stmchi+1/stpchi)/2
        gamma2=(1/stmchi-1/stpchi)/2      
        if(mi0.eq.0) then
      
         sigma=exp(-0.5*(x/eta1)**2)/eta1/sqrt2pi 
         if(thb.lt.cr.chi) then  
           MOSREF=(1+cgeom)/2*
     1       (1-exp(-qkin*cr.dw*gamma1*cr.frame.size(3)/10.*2*sigma))      ! Laue
         else  
           if (gamma2.ne.0) then                                      ! Bragg
             MOSREF=(1+cgeom)/
     1      (1+cgeom/tanh(gamma2*qkin*cr.dw*cr.frame.size(3)/10.*sigma))    
           else
             MOSREF=
     1       1/(1+1/(gamma1*qkin*cr.dw*cr.frame.size(3)/10.*sigma))
           endif
         endif 
        
        else   ! absorbing case
          sigma=qkin*cr.dw/mi0*exp(-0.5*(x/eta1)**2)/eta1/sqrt2pi    
          rb=sqrt(1+2*sigma+(sigma*cgeom)**2)
          rl=sqrt((1+2*sigma)*cgeom**2+sigma**2)
          if(thb.lt.cr.chi) then  
           dum=sigma*(1+cgeom)*exp(-(1+sigma)*cr.frame.size(3)/
     1        10.*mi0*gamma1)
           MOSREF=dum*sinh(cr.frame.size(3)/10*mi0*gamma1*rl)/rl        ! Laue
          else 
           MOSREF= sigma*(1+cgeom)/(1+sigma+rb/tanh(cr.frame.size(3)/
     1      10.*mi0*gamma1*rb))   ! Bragg
          endif
        endif 
      endif               
      end  
        
#     ---------------------------------------------------
      SUBROUTINE MAKEUPCASE(line)
#     converts LINE to uppercase     
#     ---------------------------------------------------
      implicit none
      character*(*) line
      integer*4 l,i
      l=len(line)
      do i=1,l
         if((line(i:i).ge. 'a').and.(line(i:i).le. 'z')) then
            line(i:i)=char(ichar(line(i:i))-32)
         endif  
      end do 
      end