src/ness/ness_mcryst.f

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

Source module last modified on Wed, 30 Mar 2005, 0:36;
HTML image of Fortran source automatically generated by for2html on Mon, 23 May 2005, 21:29.


#//////////////////////////////////////////////////////////////////////
#////                                                              //// 
#////  NEutron Scattering Simulation - v.2.0, (c) J.Saroun, 2000   ////
#////                                                              //// 
#//////////////////////////////////////////////////////////////////////
#////
#////  Subroutines describing objects - CRYSTAL ARRAY
#////  
#////                          
#//////////////////////////////////////////////////////////////////////
# ---------------------------------------
      SUBROUTINE CRYST_INIT2(cr)
# ---------------------------------------
      implicit none
      INCLUDE 'const.inc'
      INCLUDE 'structures.inc'                        
      record /CRYSTAL/ cr
      real*8 z,GETQKIN,GETREFDYN,GETMI
      integer*4 i,j

#///  THB,DHKL,etc.. must be specified before !

      call SLIT_INIT(cr.frame)
        
      cr.gtot=2*pi/abs(cr.dhkl)
      cr.g(1)=cr.gtot*sin(cr.chi)
      cr.g(2)=0.
      cr.g(3)=cr.gtot*cos(cr.chi)
      cr.stmch=sin(cr.thb-cr.chi)
      cr.ctmch=cos(cr.thb-cr.chi)
      cr.lambda=2.*cr.dhkl*sin(cr.thb)
      do i=1,3
        cr.mapg(i)=.false.
        do j=1,3 
          cr.dg_dr(i,j)=0
        enddo
      enddo               

# G-gradient for elastically bent crystal      
      if(cr.hmos.lt.sec) then 
          cr.dg_dr(1,1)=-cos(cr.chi)*cr.gtot*cr.rh
          cr.dg_dr(1,3)=sin(cr.chi)*cr.gtot*cr.rh
          cr.dg_dr(3,1)=sin(cr.chi)*cr.gtot*cr.rh
          cr.dg_dr(2,2)=0     ! no vertical bending
          cr.dg_dr(3,3)=-cr.poi*cos(cr.chi)*cr.gtot*cr.rh
          cr.mapg(1)=.true.        
          cr.mapg(3)=.true.        
      endif
# d-gradient
      if(cr.dgr.ne.0.) then         
          z=1.d-4*cr.gtot*cr.dgr
          cr.dg_dr(1,1)=cr.dg_dr(1,1)+z*cos(cr.dga+cr.chi)
          cr.dg_dr(1,3)=cr.dg_dr(1,3)-z*sin(cr.dga+cr.chi)
          cr.dg_dr(3,1)=cr.dg_dr(3,1)-z*sin(cr.dga+cr.chi)
          cr.dg_dr(3,3)=cr.dg_dr(3,3)-z*cos(cr.dga+cr.chi)
          cr.mapg(1)=.true.        
          cr.mapg(3)=.true.        
      endif                  
# unit vector |- to G
      cr.gama(1)=cos(cr.chi)
      cr.gama(3)=-sin(cr.chi)
      cr.gama(2)=0 

      cr.qhkl=GETQKIN(cr,cr.lambda)   ! kin. reflectivity
      cr.dext=cr.dhkl/cr.lambda*sqrt(4*pi/cr.qml)  ! extinction length
      z=cr.dlam/cr.dext
      if(z.gt.1d-5.and.cr.hmos.gt.sec) then 
         cr.ext1=tanh(z)/z          ! primary extinction
      else
         cr.ext1=1.
      endif         
      cr.mi=GETMI(cr,cr.lambda,300.d0)   ! absorption coefficient
      cr.ref=GETREFDYN(cr,cr.lambda)     ! peak reflectivity (no mosaic)         
      cr.delta=cr.qhkl*cr.dext*1d-4/pi    ! Darwin box width
      if(cr.hmos.gt.sec) cr.hmos=max(cr.hmos,cr.delta)
      
      end


#        -------------------------------------------
      logical FUNCTION CRYST_GO2(cryst,neui,neuf)
#        -------------------------------------------
      implicit none
      INCLUDE 'const.inc'
      INCLUDE 'structures.inc'
      INCLUDE 'ness_common.inc'
      INCLUDE 'randvars.inc'
                        
      record /CRYSTAL/ cryst
      record /NEUTRON/ neui,neuf
      real*8 v(3),k(3),r(3),kag(3)
      real*8 pp,dt,h,tin,tout
      real*8 CRYST_ARRAY
      integer*4 i

#      logical*4 ref  
#      ref=(cryst.frame.size(2).eq.80.and.cryst.frame.count.gt.1000)      
1     format(a10,2x,i8,2x,10(1x,g12.6))
#2     format(a10,(1x,g12.6),5x,$)
      if(dbgref) then
#        write(*,1) 'XRND ',(XRND(i),i=1,5)
      write(*,1) cryst.frame.name,cryst.frame.count,neui.r,neui.k,neui.p
      endif  

#/// if nh=0 then accept neutrons without transformations
      if (cryst.nh.eq.0) then           
        neuf=neui
        dt=(cryst.frame.dist-neui.r(3))/neui.k(3)
        do i=1,3
              neuf.r(i)=neui.r(i)+dt*neui.k(i)
        end do
        neuf.t=neui.t+dt/hovm
        goto 10
      endif   

      if (cryst.mag*neui.s.lt.0) goto 99      

#      if (ref) then
#         call wrtneu(NEUI)
#         write(*,1) 'dAlpha1: ',atan(NEUI.K(1)/NEUI.K(3))/deg,NEUI.P
#      endif   
        
      call SLIT_PRE(cryst.frame,neui.r,neui.k,v,k)
      v(1)=v(1)+cmx              
# test right barrier
      if (cbar.gt.0) then
         h=v(3)-(cryst.frame.size(1)/2.+v(1)+40.)/k(1)*k(3)
         if (abs(h).lt.cbar.and.k(1).gt.0) goto 99
      endif
      
#/// move to the entry point          
      call CR_BORDER(cryst,v,k,tin,tout)            
      if (tin.ge.tout) goto 99      ! No intersection with the crystal
      do i=1,3
         v(i)=v(i)+(tin+1.d-7*(tout-tin))*k(i)
      end do

#      if (ref) write(*,1) 'start: ',(V(I),i=1,3),NEUI.P

      pp=CRYST_ARRAY(cryst,v,k,r,kag,dt)
      
      if (pp.gt.1.d-4) then
# test right barrier
          if (cbar.gt.0) then
            h=r(3)-(cryst.frame.size(1)/2.+r(1)+40.)/kag(1)*kag(3)
            if (abs(h).lt.cbar.and.kag(1).lt.0) goto 99
          endif
# transform to local axis coordinate and return
          r(1)=r(1)-cmx                                 
          neuf.s=neui.s
          neuf.t=neuf.t+dt/hovm
          neuf.p=neui.p*pp
#      if (ref) write(*,1) 'end: ',(R(I),i=1,3),NEUF.P
          call SLIT_POST(cryst.frame,r,kag,neuf.r,neuf.k)
#      if (ref) then
#         call wrtneu(NEUF)
#         write(*,1) 'dAlpha2: ',
#     *              (atan(NEUF.K(1)/NEUF.K(3))+2*CRYST.THB)/deg,NEUF.P
#      endif   
          cryst.frame.count=cryst.frame.count+1
          goto 10
      else
          goto 99
      endif
        
10    CRYST_GO2=.true.
#      if (dbgref)  then
#      write(*,1) CRYST.FRAME.NAME,CRYST.FRAME.COUNT,NEUF.R,NEUF.K,NEUF.P
#      endif
       return
      
99    neuf.p=0
#      if (dbgref)  then
#      write(*,1) CRYST.FRAME.NAME,CRYST.FRAME.COUNT,NEUF.P
#      endif
      CRYST_GO2=.false.
      return        
      end        
        

#---------------------------------------------------------------
      real*8 FUNCTION HMOS_DIST(x)
# Distribution of horizontal angular deviation of mosaic blocks
# --------------------------------------------------------------
      implicit none
      INCLUDE 'const.inc'
      INCLUDE 'structures.inc'
      INCLUDE 'ness_common.inc'
      real*8 x,sqrt2pi,sqrt2ln2,dd0
      parameter (sqrt2pi=2.506628274631, sqrt2ln2=1.177410023, 
     *           dd0=0.270347526)      

      if (mdist.eq.1) then        
# pseudo-Voigt
        HMOS_DIST=0.5*dd0/(1. + (x/sqrt2ln2)**2)+
     *   0.5*exp(-0.5*x**2)/sqrt2pi
      else if (mdist.eq.2) then        
# Lorenz
        HMOS_DIST=dd0/(1. + (x/sqrt2ln2)**2)
      else if (mdist.eq.3) then        
# Rectangle
        if (abs(x).le.0.5) then
          HMOS_DIST=1.
        else
          HMOS_DIST=0. 
        endif   
      else 
# Gauss
          HMOS_DIST=exp(-0.5*x**2)/sqrt2pi
      endif
      end 

#-----------------------------------------------------------
      real*8 FUNCTION CRYST_ARRAY(cr,r1,k1,r2,k2,dt)
# Transmission function for mosaic crystal 
# (eta>")  
# simplified version (Darwin width=0)    
#-----------------------------------------------------------
      implicit none
      INCLUDE 'const.inc'
      INCLUDE 'structures.inc'
      INCLUDE 'crystal.inc'
      real*8 eps
      parameter (eps=1d-5)
      record /CRYSTAL/ cr
      real*8 r1(3),k1(3),r2(3),k2(3),dt
      real*8 g(3)
      integer*4 i
      real*8 pp,tin,tout
      real*8 kk,k0,qhkl,z,mi

#      logical*4 ref      
#      ref=(cr.frame.size(2).eq.81.and.cr.frame.count.gt.1000)      
1     format(a,10(1x,g12.6))
      

      pp=1.d0
#/// calculate values, which are constant in the course of the iteration      
      kk=k1(1)**2+k1(2)**2+k1(3)**2
      k0=sqrt(kk)      
      qhkl=cr.qhkl/10.*cr.dw*cr.ext1
      mi=cr.mi/10.
      path=0.d0
      tof=0.d0
      call CR_BORDER(cr,r1,k1,tin,tout)
            
      if (tin.ge.tout) goto 99      ! No intersection with the crystal
#/// move to the entry point = 1st segment sector entry         
      do i=1,3
         r2(i)=r1(i)+(tin+1.d-7*(tout-tin))*k1(i)
      end do
      
# ****  Begin of multiple scattering iteration cycle  ****

#// generate random walk step in K direction
#---------------------------------------------------
50    call WALKSTEP(cr,1,r2,k1,k0,qhkl,g,pp,tout)     
      if (pp.lt.eps) goto 99    ! don't follow trajectories with low probability         
#// move to the new reflection point and set K=K+G
      do i=1,3
         k2(i)=k1(i)+g(i)
      end do
      
      z=sqrt(kk/(k2(1)**2+k2(2)**2+k2(3)**2))   ! set |K2| = |K1|
      do i=1,3
           k2(i)=k2(i)*z
      end do
      
#      if (dbgref) then
#        write(*,1) 'dTheta1: ',atan(G(1)/G(3))/deg
#        write(*,1) 'step1:  ',(r2(i),i=1,3),PP,
#     *              (atan(-k1(3)/k1(1))-CR.THB)/deg,
#     *              (atan(k2(3)/k2(1))-CR.THB)/deg
#        pause
#      endif  

#// generate 2nd random walk step in K+G direction
#-------------------------------------------------
      call WALKSTEP(cr,-1,r2,k2,k0,qhkl,g,pp,tout)      
      
      if (pp.lt.eps) goto 99
      if (tout.le.0) then    
        pp=pp*exp(-mi*path)
      else       
        do i=1,3
           k1(i)=k2(i)+g(i)
        end do
        z=sqrt(kk/(k1(1)**2+k1(2)**2+k1(3)**2))   ! set |K2| = |K1|
        do i=1,3
           k1(i)=k1(i)*z
        end do        
        
#        if (dbgref) then
#          write(*,1) 'dTheta2: ',atan(G(1)/G(3))/deg
#          write(*,1) 'step2:  ',(r2(i),i=1,3),PP,
#     *              (atan(-k1(3)/k1(1))-CR.THB)/deg,
#     *              (atan(k2(3)/k2(1))-CR.THB)/deg
#          pause
#        endif
        
         goto 50   ! repeat cycle
      endif
        
# ****  End of multiple scattering cycle  ****
      
#      if (ref) then
#        write(*,1) 'end:  ',(r2(i),i=1,3),pp
#        write(*,1) 'path: ',path,tof*k0,pp
#        pause  
#      endif
      
      dt=tof
      CRYST_ARRAY=pp
      return
99    CRYST_ARRAY=0.d0
#      IF (DBGREF) THEN
#        WRITE(*,1) 'EX:   ',PP,TOUT*K1(3)
#        pause     
#      ENDIF
      end         


#--------------------------------------------------------------
      SUBROUTINE WALKSTEP(cr,dir,r,k,k0,q,g,pp,tout)
# Generate random walk step in the crystal
#//   Q = Qhkl*DW*Ext1
#//   K0=|K|
#//   R,K position and K-vector
#//   DIR=1,-1 for directions K,K+G
#//   TOUT  .Time to reach assembly exit
#//   G(3) ...Local diffraction vector at R
#//   PP ... weight to be attributed to the event
#--------------------------------------------------------------
      implicit none
      INCLUDE 'structures.inc'
      INCLUDE 'crystal.inc'
      record /CRYSTAL/ cr
      real*8 eps
      parameter(eps=1d-5)
      integer*4 i,m,dir,i0(3)
      real*8 z,phi,p,dt,k0,q,tout,pp
      real*8 r(3),k(3),r0(3),g(3)
      real*4 RAN1,GASDEV1
      
#      logical*4 ref      
#      ref=(cr.frame.size(2).eq.81.and.cr.frame.count.gt.1000)      
#1     format(a,10(1x,g12.6))


      
# add random vertical mosaic angle 
      phi=cr.vmos*GASDEV1(1,0.,3.)       
# trace through segments and find times, scatt. probabilities, angular dev. etc..
      call CR_BORDERS(cr,dir,r,k,k0,q,phi) 
      
      if (nseg.le.0) goto 99
            
      tout=tseg(nseg)  ! time to the last segment exit
      p=pseg(nseg)    ! scattering probability up to the assembly exit 
      if(dir.lt.0) p=1.  ! for K+G direction, consider the possibility of leaving the assembly  
      if (p*pp.lt.eps) goto 99  ! don't follow trajectories with low probability
# Select a segment for next reflection according to the scattering probabilities
      z=p*RAN1()
      m=1
      do while (m.le.nseg.and.z.ge.pseg(m))
          m=m+1
      end do
      
# No reflection in any segment
      if(m.gt.nseg) then
        if(dir.lt.0) then
            goto 90  ! go to the last segment exit and return
        else    
            goto 99  ! event stopped (no reflection) 
        endif   
      endif          

# Find point of reflection inside the M-th segment 
      call SEGTRACE(cr,k0,q,m,dt)  
      
      if (dt.le.0) goto 99

# accumulate flight-path through the material       
      do i=1,m-1
         path=path+(tseg(i)-tseg0(i))*k0
      enddo
      path=path+dt*k0

# accumulate time-of-flight     
      tof=tof+dt+tseg0(m)
                  
# Move to the point of reflection
      do i=1,3
         r(i)=r(i)+(dt+tseg0(m))*k(i)
      enddo      
      
#      IF(REF) THEN       
#       write(*,1) 'DIR: ',DIR
#       write(*,1) 'DT : ',DT*K0,TSEG0(M)*K0       
#       write(*,1) 'R1: ',(R(I)-(DT+TSEG0(M))*K(I),I=1,3)
#       write(*,1) 'R2: ',(R(I),I=1,3)
#      pause
#      ENDIF
      

# Coordinates of the segment centre
      call SEGCOORD(cr,r,r0,i0)
        
#      if (ref)  write(*,1) 'segment  : ',I0(1)
# Calculate local G-vector 
      call LOCALG(cr,dir,r,r0,i0,alpha(m)+grad(m)*dt,phi,g)
      
#      IF (REF) THEN        
#        Z=v3xv3(G,G)
#        DO I=1,3
#          AT(I)=K(I)+G(I)
#        ENDDO
#        write(*,1) 'K  : ',(K(I),I=1,3)
#        write(*,1) 'G  : ',(G(I),I=1,3)
#        write(*,1) 'DELTA: ',2.*v3xv3(AT,G)-Z        
#        CALL LOCALG(CR,-DIR,R,R0,I0,ALPHA(M)+GRAD(M)*DT,PHI,G0)
#        write(*,1) 'K  : ',(AT(I),I=1,3)
#        write(*,1) 'G  : ',(G0(I),I=1,3)
#        Z=v3xv3(G0,G0)
#        DO I=1,3
#          AT(I)=AT(I)+G0(I)
#        ENDDO
#        write(*,1) 'DELTA: ',2.*v3xv3(AT,G0)-Z
#      ENDIF
      
      pp=pp*p
      return
            
# No reflection when going along K+G (DIR<0) 
90    do i=1,3
         r(i)=r(i)+tseg(nseg)*k(i)
      enddo
      tof=tof+tseg(nseg)
      do i=1,nseg
          path=path+(tseg(i)-tseg0(i))*k0
      enddo
      tout=0.d0
      
      return
      
# Left without reflection 
99    pp=0.d0
      tout=0.d0
            
      end

# ---------------------------------------------------------------------C
      SUBROUTINE CR_BORDERS(cr,dir,r,k,k0,q,phi)
# Traces through all segments along the neutron path
# PHI= vertical mosaic tilt. angle 
# Q = Qhkl*DW*Ext1
# K0=|K|
# R,K position and K-vector
# DIR=1,-1 for directions K,K+G
# IN /CRBORDERS/ returns for each segment crossed:
# TSEG    time to the I-th segment exit, started at the assembly entry
# TSEG0   time to the I-th segment entry, -"-
# ALPHA   dThetaB, deviation from Bragg angle 
# GRAD    grad(ThetaB), gradient of Bragg angle along K 
# PSEG    scattering probability
# NSEG    number of crossed segments
# PATH    accumulates path-length through the material
# ---------------------------------------------------------------------C
      implicit none
      INCLUDE 'structures.inc'
      INCLUDE 'crystal.inc'
      record /CRYSTAL/ cr
      integer*4 i,dir,j,i0(3)
      real*8 kk,k0,t,dt,q,phi,tin1,tin2,tout1,tout2     
      real*8 r(3),k(3),v(3),v0(3)
      real*8 limit
      parameter(limit=0.5d0)
#      LOGICAL*4 REF
#      REF=(CR.FRAME.COUNT.GE.10000)      
      
#1     FORMAT(a11,6(1x,G16.10))

      j=0
      tseg(0)=0.d0
      tseg0(0)=0.d0
      pseg(0)=0.d0
      t=0.d0  ! T measures time-fo-flight
      kk=k0**2
      do i=1,3
         v(i)=r(i)
      end do
#      if(ref)  write(*,1) 'START AT: ',(V(I),I=1,3),DIR
      
      do while ((abs(v(1)/cr.frame.size(1)).lt.limit).and.
     *          (abs(v(2)/cr.frame.size(2)).lt.limit).and.
     *          (abs(v(3)/cr.frame.size(3)).lt.limit).and.j.lt.mseg)
      
        call SEGCOORD(cr,v,v0,i0)  ! get segment coordinates
        
#       if(ref)  write(*,1) 'CENTER: ',(V0(I),I=1,3),(I0(I),I=1,3)
         
# Get entry (TIN) and exit (TOUT) times of a neutron moving along K and starting at V
# (1) .. crossing segment sector
# (2) .. crossing the segment itself (differs by gaps between the segments, CR.DH ...)
        call SEGCROSS(cr,v,k,v0,tin1,tin2,tout1,tout2)

#       if (ref) then
#         write(*,1) 'CROSS IN : ',(TIN1*K(I),I=1,3),(TIN2*K(I),I=1,3)
#         write(*,1) 'CROSS OUT: ',(TOUT1*K(I),I=1,3),(TOUT2*K(I),I=1,3)
#       endif  
                          
# Count only intersected segments
        if((tout2-tin2)*k0.gt.1.d-3.and.tout2.gt.1.d-10) then  ! must start before or inside the segment
          j=j+1        
# Move to the segment entry, if not already inside   
          if (tin2.gt.0.) then 
            do i=1,3
              v(i)=v(i)+tin2*k(i)
            enddo 
            t=t+tin2
          else
            tin2=0.d0
          endif
          tseg0(j)=t      ! J-th segment starting time (with respect to the orig. position R)    
          tseg(j)=t+tout2-tin2   ! J-th segment exit time     
                          
# Calculate total scattering probability for the segment, starting at V          
          call SEGSCATT(cr,j,dir,v,k,phi,q,tseg(j)-tseg0(j))
          
#       if(ref) THEN
#          write(*,1) 'SEG POS : ',(V(I),I=1,3)
#          write(*,1) 'SEG EXIT: ',(V(I)+K(I)*(TSEG(J)-TSEG0(J)),I=1,3)
#          write(*,1) 'ALPHA: ',ALPHA(J),GRAD(J),-ALPHA(J)/GRAD(J)*K0
#          write(*,1) 'PSEG : ',PSEG0(J),PSEG(J)
#       endif

        else
          tin2=0.d0        
        endif

# Move to the entry of the next segment sector        
        dt=1.000001d0*(tout1-tin2)   ! move slightly behind to have it realy INSIDE the next sector
        do i=1,3
           v(i)=v(i)+dt*k(i)
        enddo         
        t=t+dt   
           
#        if(ref) write(*,1) 'NEXT SEG: ',J,(V(I),I=1,3) 
      
      enddo
                  
      nseg=j
      
      end
      

# -------------------------------------------------------------
      SUBROUTINE LOCALG(cr,dir,r,r0,i0,eta,phi,g)
# Calculate local G-vector at R
# I0(3) segment coordinates
# R0(3) segment center physical coordinates
# ETA .. horizontal tilt angle of mosaic domain
# PHI .. vertical tilt angle of mosaic domain
# -------------------------------------------------------------
      implicit none
      INCLUDE 'structures.inc'
      record /CRYSTAL/ cr
      integer*4 i,dir,i1,i0(3)
      real*8 r(3),r0(3),g(3),w(3),at(3),g0(3)
      real*8 z,eta,phi,gabs
      
# Calculate local G-vector (only G-gradient) 
      do i=1,3
          g0(i)=cr.g(i) 
          if (cr.mapg(i)) then
              w(i)=0.d0
            do i1=1,3
                w(i)=w(i)+cr.dg_dr(i,i1)*(r(i1)-r0(i1))
            enddo
            g0(i)=g0(i)+w(i)
          endif  
      end do
      gabs=sqrt(g0(1)**2+g0(2)**2+g0(3)**2)

# Add segment tilt angle and vertical mosaic spread 
      call SEGTILT(cr,i0,at)
      g(1)=g0(1)-g0(3)*(at(1)+at(3))  
      g(3)=g0(3)+g0(1)*(at(1)+at(3)) 
      g(2)=g0(2)-g0(3)*at(2)+gabs*phi

# Add the angle of the mosaic block
      do i=1,3
         g(i)=g(i)+gabs*cr.gama(i)*eta
      end do
# Renormalize
      z=gabs/sqrt(g(1)**2+g(2)**2+g(3)**2)
      do i=1,3
         g(i)=dir*g(i)*z
      enddo
      end


# -------------------------------------------------------------
      SUBROUTINE CR_BORDER(cr,r,k,tin,tout)
#     Returns times of intersection with the crystal assembly borders, 
#     started at current position R and measured along K.
#     All in crystal local coordinate.
#  !! Time is in units [sec*h/m] everywhere   i.e. length=time*K
#--------------------------------------------------------------
      implicit none
      INCLUDE 'structures.inc'
#      INCLUDE 'crystal.inc'
      record /CRYSTAL/ cr
        
      real*8 r(3),k(3),tin,tout,t1(3),t2(3),dum
      integer*4 i
        
      do i=1,3
      if (abs(k(i)).gt.1.0d-8) then
               t2(i)=(cr.frame.size(i)/2.d0 - r(i))/k(i)
        t1(i)=(-cr.frame.size(i)/2.d0 - r(i))/k(i)
        if (t1(i).gt.t2(i)) then
          dum=t1(i)
          t1(i)=t2(i)
          t2(i)=dum
        endif            
      else
          t2(i)=1.0d30
          t1(i)=-1.0d30
      endif
      end do
      tin=max(t1(1),t1(2),t1(3))
      tout=min(t2(1),t2(2),t2(3))
      if (tin.gt.tout) then
         tin=1d30
         tout=1d30
      endif
      end            



# ----------------------------------------------------------
      SUBROUTINE SEGTILT(cr,i0,at)
#     return tilt angles
# ----------------------------------------------------------
      implicit none
      INCLUDE 'structures.inc'
      record /CRYSTAL/ cr
      real*8 at(3),da
      integer*4 i0(3)
      
      if(i0(1).gt.0.and.i0(1).le.cr.nh) then
           da=cr.frame.size(1)*cr.rh/cr.nh        
           at(1)=(-(cr.nh-1)/2.d0+i0(1)-1.d0)*da
      else 
           at(1)=0.d0
      endif
      if(i0(2).gt.0.and.i0(2).le.cr.nv) then
           da=cr.frame.size(2)*cr.rv/cr.nv        
           at(2)=(-(cr.nv-1)/2.d0+i0(2)-1.d0)*da
      else 
           at(2)=0.d0
      endif
      if(i0(3).gt.0.and.i0(3).le.cr.nb) then
           da=cr.frame.size(3)*cr.rb/cr.nb        
           at(3)=(-(cr.nb-1)/2.d0+i0(3)-1.d0)*da
      else 
           at(3)=0.d0
      endif
      
#      IF(CR.HMOS.GT.0) THEN
#      write(*,*) CR.FRAME.NAME,CR.NH,CR.NV,CR.NB
#1     format (3(G10.4,2x))      
#      write(*,*) (I0(I),I=1,3)
#      write(*,1) (AT(I),I=1,3)
#     pause      
#      ENDIF
      end
      
      
# ---------------------------------------------------------------------------
      SUBROUTINE SEGCOORD(cr,r,r0,i0)
#     return coordinates of the segment R0, in which the particle at R resides
# ---------------------------------------------------------------------------
      implicit none
      INCLUDE 'structures.inc'
      record /CRYSTAL/ cr
      real*8 r(3),r0(3),half
      integer*4 i0(3) 
      parameter(half=0.5d0)

# ih,iv,ib are the integre-coordinates of the closest segment
      i0(1)=nint((r(1)/cr.frame.size(1)+half)*cr.nh+half)
      i0(2)=nint((r(2)/cr.frame.size(2)+half)*cr.nv+half)
      i0(3)=nint((r(3)/cr.frame.size(3)+half)*cr.nb+half)
# get physical coordinates of the segment center
      r0(1)=cr.frame.size(1)*(1.d0*(i0(1)-half)/cr.nh-half)
      r0(2)=cr.frame.size(2)*(1.d0*(i0(2)-half)/cr.nv-half)
      r0(3)=cr.frame.size(3)*(1.d0*(i0(3)-half)/cr.nb-half)
      end


# ---------------------------------------------------------------------------C
      SUBROUTINE SEGCROSS(cr,r,k,r0,tin1,tin2,tout1,tout2)
#     return times when the particle (R,K) crossess the boarders
# of the segment sector (TIN1,TOUT1) and the segment itself (TIN2,TOUT2) 
# ---------------------------------------------------------------------------C
      implicit none
      INCLUDE 'structures.inc'
      record /CRYSTAL/ cr
      real*8 r(3),r0(3),k(3),tin1,tin2,tout1,tout2
      real*8 r1(3),bsz(3),t1(3),t2(3),dum
      integer*4 i 
1     format(a10,6(2x,g12.6))

# sector size
      bsz(1)=cr.frame.size(1)/cr.nh
      bsz(2)=cr.frame.size(2)/cr.nv
      bsz(3)=cr.frame.size(3)/cr.nb
# coordinates relative to segment centre
      do i=1,3      
        r1(i)=r(i)-r0(i)
      enddo
# Get entry and exit times for segment sector      
      do i=1,3      
      if (abs(k(i)).gt.1.0d-10) then
               t2(i)=(bsz(i)/2.-r1(i))/k(i)
        t1(i)=(-bsz(i)/2.-r1(i))/k(i)
        if (t1(i).gt.t2(i)) then
          dum=t1(i)
          t1(i)=t2(i)
          t2(i)=dum
        endif            
      else
          t2(i)=1.0d30
          t1(i)=-1.0d30
      endif
      enddo
      tin1=max(t1(1),t1(2),t1(3))
      tout1=min(t2(1),t2(2),t2(3))
# segment size
      bsz(1)=cr.frame.size(1)/cr.nh-cr.dh
      bsz(2)=cr.frame.size(2)/cr.nv-cr.dv
      bsz(3)=cr.frame.size(3)/cr.nb-cr.db
# Get entry and exit times for the segment itself      
      do i=1,3      
      if (abs(k(i)).gt.1.0d-10) then
               t2(i)=(bsz(i)/2.d0-r1(i))/k(i)
        t1(i)=(-bsz(i)/2.d0-r1(i))/k(i)
        if (t1(i).gt.t2(i)) then
          dum=t1(i)
          t1(i)=t2(i)
          t2(i)=dum
        endif            
      else
          t2(i)=1.0d30
          t1(i)=-1.0d30
      endif
      enddo

#       if (R(1).LT.-8.5) then
#         write(*,1) 'R1: ',(R1(I),I=1,3),(BSZ(I),I=1,3)
#         write(*,1) 'T1: ',(T1(I),I=1,3)
#         write(*,1) 'T2: ',(T2(I),I=1,3)
#         pause
#       endif

      tin2=max(t1(1),t1(2),t1(3))
      tout2=min(t2(1),t2(2),t2(3))

      end

#-----------------------------------------------------------------------------
      SUBROUTINE SEGSCATT(cr,j,dir,r,k,phi,q,dt)
# Calculate scattering probability for J-th segment on the path along K(3)
# PHI is the vertical mosaic angle
# Q .. kinematical reflectivity (incl. DW) 
# R is the coordinate of the neutron in the assembly
# DT ... time-of-flight along K(3) through the segment 
#-----------------------------------------------------------------------------
      implicit none
      INCLUDE 'structures.inc'
      INCLUDE 'crystal.inc'
      record /CRYSTAL/ cr
      real*8 sec
      parameter (sec=4.85d-6)
      integer*4 dir,j,i,i1,i0(3),l
      real*8 dt,a,b,z,z1,z2,sigma,kk,k0,gabs,q,phi,dalpha
      real*8 kag(3),k(3),g(3),r(3),r0(3),v(3)
      real*8 HMOS_DIST,ERF
#1     FORMAT(a10,10(1x,G12.6))
      
      if (dt.le.0) then
          pseg0(j)=0.d0
          pseg(j)=pseg(j-1)
          return
      endif      
      
      kk=0.
      do i=1,3
         kk=kk+k(i)**2
      end do
      k0=sqrt(kk)
      
# get coordinates of the segment centre
      call SEGCOORD(cr,r,r0,i0)
      
      call LOCALG(cr,dir,r,r0,i0,0.d0,phi,g)
      
# get angular deviation from the Bragg condition (=alpha)            
      gabs=0.     
      do i=1,3
         kag(i)=k(i)+g(i)
         gabs=gabs+g(i)**2
      end do
      gabs=sqrt(gabs)
      a=kag(1)**2+ kag(2)**2+kag(3)**2-kk     ! (K+G)^2-K^2 
      b=dir*gabs*(kag(1)*cr.gama(1)+kag(2)*cr.gama(2)+
     *       kag(3)*cr.gama(3))   

      grad(j)=0.d0
      do i=1,3
        if (cr.mapg(i)) then
          z=0.d0
          do i1=1,3
            z=z+dir*cr.dg_dr(i,i1)*k(i1)
          enddo
          grad(j)=grad(j)+kag(i)*z
        endif  
      end do
# get gradient of angular deviation ALPHA            
      grad(j)=-grad(j)/b
# get angular deviation from the Bragg condition (=alpha)            
      alpha(j)=-a/(2.*b)
      
#      write(*,1) 'ALPHA: ',alpha(J)
# 2nd order correction for large bent crystals
      if(abs(grad(j)).ge.1.d-6*k0*q.and.cr.hmos.le.sec) then
      dalpha=alpha(j)
      l=10  ! max. L iterations
      do while (abs(dalpha).gt.1d-7.and.l.gt.0)
        l=l-1
        do i=1,3
          v(i)=r(i)-alpha(j)/grad(j)*k(i)
        enddo
        call LOCALG(cr,dir,v,r0,i0,0.d0,phi,g)
        gabs=0.     
        do i=1,3
          kag(i)=k(i)+g(i)
          gabs=gabs+g(i)**2
        end do
        gabs=sqrt(gabs)
        a=kag(1)**2+ kag(2)**2+kag(3)**2-kk     ! (K+G)^2-K^2 
        dalpha=-a/(2.*b)
        alpha(j)=alpha(j)+dalpha 
#        write(*,1) 'DALPHA: ',-a/(2.*b)
      enddo
      endif
#      pause


#// add random angle within Darwin box width
#        alpha(J)=alpha(J)+(RAN1(1)-0.5)*Q*CR.DEXT*1D-3/PI

# MOSAIC CRYSTAL
      if(abs(grad(j)).lt.1.d-6*k0*q.and.cr.hmos.gt.sec) then
        sigma=k0*q*HMOS_DIST(alpha(j)/cr.hmos)/cr.hmos*dt
      else  
# MOSAIC CRYSTAL WITH G-GRADIENT
        if(cr.hmos.gt.sec) then
           z1=ERF((alpha(j)+grad(j)*dt)/cr.hmos,0)
           z2=ERF(alpha(j)/cr.hmos,0)
           sigma=k0*q*(z1-z2)/grad(j)
           if(j.eq.1.and.abs(alpha(j)).lt.sec) sigma=0.d0  ! no back refl. in the same segment
# ONLY G-GRADIENT
        else   
           z1=-alpha(j)/grad(j)
           if(z1.gt.0.and.z1.lt.dt) then  
             sigma=k0*q/abs(grad(j))
             if(j.eq.1.and.abs(alpha(j)).lt.sec) sigma=0.d0  ! no back refl. in the same segment   
           else
             sigma=0.d0
           endif
        endif                              
      endif        
      if(sigma.gt.14) then
          pseg0(j)=1.d0
          pseg(j)=1.d0
      else  
            pseg0(j)=1.d0-exp(-sigma)
          pseg(j)=1.d0-(1.d0-pseg(j-1))*(1.d0-pseg0(j))
      endif  
      
      end
      

#-------------------------------------------------------------------------
      SUBROUTINE SEGTRACE(cr,k0,q,m,dt)
# return random step size to a point of reflection inside the M-th segment
#-------------------------------------------------------------------------
      implicit none
      INCLUDE 'structures.inc'
      INCLUDE 'crystal.inc'
      record /CRYSTAL/ cr
      integer*4 m
      real*8 sec,ksi,z,z1,sigma,p0,dt,aa
      parameter(sec=4.85d-6)
      real*8 k0,q
      real*4 RAN1
      real*8 ERF,HMOS_DIST
     
      dt=0.d0
      
# MOSAIC CRYSTAL
      if(abs(grad(m)).lt.1.d-6*k0*q) then
          ksi=RAN1()
          p0=pseg0(m)
          aa=alpha(m)/cr.hmos
          sigma=k0*q*HMOS_DIST(aa)/cr.hmos
          dt=-log(1-ksi*p0)/sigma
      else
# MOSAIC CRYSTAL WITH D-GRADIENT
        if(cr.hmos.gt.sec) then
          ksi=RAN1()
          p0=pseg0(m)
          z=grad(m)*log(1-ksi*p0)/(k0*q)
          aa=alpha(m)/cr.hmos
#  ERF shoud be used only on (-inf.;0) because of num. precision
          if(aa.gt.0) then
            z1=ERF(-aa,0)
            if(1.-z1-z.gt.0.5d0) then
              dt=-ERF(z1+z,1)-aa
            else
              dt=ERF(1.-z1-z,1)-aa
            endif
          else
            z1=ERF(aa,0)
            if(z1-z.gt.0.5d0) then
              dt=-ERF(1.-z1+z,1)-aa
            else
              dt=ERF(z1-z,1)-aa
            endif             
          endif
           dt=dt*cr.hmos/grad(m)                     
# ONLY D-GRADIENT
        else  
          dt=-alpha(m)/grad(m)
        endif
      endif
      
      end
      
#-----------------------------------------------------
      SUBROUTINE ERF_INIT(f,amin,amax)
# Calculate lookup tables for ERF function
#-----------------------------------------------------
      implicit none
      integer*4 i,j
      real*8 sum,z1,z2,z3,a,b,det,x1i,xj,amin,amax
      integer*4 dim
      parameter(dim=1025)
      real*8 xmin,dx,y(dim),xmin1,dx1,y1(dim)
      common /erfcom/ xmin,dx,y,xmin1,dx1,y1      
      real*8 f
      external f

# Generate cumulative function from F(X)
      sum=0.
      dx=(amax-amin)/(dim-1)
      xmin=amin
      y(1)=0.
      do i=1,dim-1
        z1=xmin+(i-1)*dx
        z2=z1+dx/2
        z3=z1+dx  
        sum=sum+f(z1)+4*f(z2)+f(z3)
        y(i+1)=sum
      end do  
      do i=1,dim
        y(i)=y(i)/y(dim)
      enddo
# Generate inverse cumulative function
      dx1=1.d+0/(dim-1)
      xmin1=y(1)
      y1(1)=xmin
      y1(dim)=xmin+(dim-1)*dx
      j=1
      do i=2,dim-1
        x1i=xmin1+(i-1)*dx1
        do while (y(j).lt.x1i.and.j.lt.dim-1)
          j=j+1
        end do
10      xj=xmin+(j-1)*dx
        a=(y(j+1)+y(j-1)-2*y(j))/2
        b=(y(j+1)-y(j-1))/2
        if (abs(a).lt.1d-30) then
           j=j-1
           goto 10
        else
          det=b**2-4*a*(y(j)-x1i)
          if (det.le.0) then
            write(*,*)  'Error in ERF_INIT: ',det,a,b
            pause
          endif
          z1=xj+dx*(-b+sqrt(det))/2/a
          z2=xj+dx*(-b-sqrt(det))/2/a
          if (abs(z2-xj).lt.abs(z1-xj)) z1=z2
          y1(i)=z1
        endif
      end do
  
      end

#----------------------------------------------------------------
      real*8 FUNCTION ERF(arg,inv)
# Return cumulative function (or inverse, if INV=0)
# Uses lookup table generated by CUM_INIT
#----------------------------------------------------------------
      implicit none
      integer*4 dim,inv,j1,j2,j,i
      real*8 arg,a,b,z,xj,det,z1,z2,arg1
      parameter(dim=1025)
      real*8 xmin,dx,y(dim),xmin1,dx1,y1(dim)
      real*8 ERF_INTERP
      common /erfcom/ xmin,dx,y,xmin1,dx1,y1      

# Cumul. function
      if (inv.ne.1) then  
        ERF=ERF_INTERP(xmin,dx,y,dim,arg)
# else Iverse cumul. function
      else
        if(arg.le.xmin1) then
           ERF=y1(1)  ! left limit
        else if(arg.ge.xmin1+(dim-1)*dx1) then
           ERF=y1(dim)   ! right limit
        else
# Find J1,J2 so that Y(J1) > A >= Y(J2) and J2=J1+1
          arg1=arg
          if(arg.gt.9.d-1) arg1=1.d0-arg
          z=(arg1-xmin1)/dx1+1
          i=int(z)
          z1=(y1(i)-xmin)/dx+1
          z2=(y1(i+1)-xmin)/dx+1   
          j1=int(z1)
          j2=int(z2)+1
          do while(j2.gt.j1+1)
            j=(j2+j1)/2
            if(y(j).lt.arg1) then
               j1=j
            else
               j2=j
            endif
          enddo
# Set J so that Y(J) is close to ARG 
          j=j1
          if(arg1-y(j1).gt.y(j2)-arg1) j=j2
# but avoid J=1 or J=DIM
          if(j.lt.2) j=j+1
          if(j.gt.dim-1) j=j-1
# interpolate quadratically between Y(J-1) and Y(J+1) 
# return inverse value for Y=ARG
10        xj=xmin+(j-1)*dx
          a=(y(j+1)+y(j-1)-2*y(j))/2
          b=(y(j+1)-y(j-1))/2
          if(abs(b).lt.1.d-30) then
             ERF=xj
          else if (abs(a).lt.1d-30) then
            ERF=xj+y(j+1)
            j=j-1
            goto 10
          else
            det=b**2-4*a*(y(j)-arg1)
            det=sqrt(det)
            z1=xj+dx*(-b+det)/2/a
            z2=xj+dx*(-b-det)/2/a
            if (abs(z2-xj).lt.abs(z1-xj)) z1=z2
            if(arg.gt.9.d-1) z1=-z1
            ERF=z1
          endif
        endif            
      endif
      end 

#-------------------------------------------------------------------
      real*8 FUNCTION ERF_INTERP(xmin,dx,y,dim,a)
# Quadratic interpolation in Y array with equidistant X
# given by XMIN,DX. MIN,MAX are values to be returned outside limits
#-------------------------------------------------------------------
      implicit none
      integer*4 dim,i
      real*8 xmin,dx,a,z,xmax
      real*8 y(dim)   

      xmax=xmin+(dim-1)*dx
      if (a.le.xmin) then 
         ERF_INTERP=y(1)
      elseif (a.ge.xmax) then         
         ERF_INTERP=y(dim)
      elseif (a.le.xmin+dx) then         
         ERF_INTERP=y(1)+(y(2)-y(2))*(a-xmin)/dx
      elseif (a.ge.xmax-dx) then         
         ERF_INTERP=y(dim-1)+(y(dim)-y(dim-1))*(a-xmax+dx)/dx
      else
        z=(a-xmin)/dx+1
        i=nint(z)
        if (i.eq.dim-1) i=i-1
        if (i.eq.2) i=3
        ERF_INTERP=y(i)+(y(i+1)-y(i-1))/2*(z-i)+
     *             (y(i+1)+y(i-1)-2*y(i))/2*(z-i)**2
      endif
      end

      
#        --------------------------------------------------
        logical*4 FUNCTION CR_INSIDE(cryst,r)
#       INSIDE function for CRSYTAL object ... takes into
#       account curved surface of a bent crystal plate.        
#     NOT USED IN THE CURRENT VERSION !
#        --------------------------------------------------        
        implicit none

        INCLUDE 'structures.inc'
      
        record /CRYSTAL/ cryst
        real*8 r(3),r0(3)
        logical*4 INSIDE
        
        
        r0(3)=r(3)-r(1)**2*cryst.rh-r(2)**2*cryst.rv
        r0(1)=r(1)
        r0(2)=r(2)
        CR_INSIDE=INSIDE(cryst.frame,r0)
        
        return
        end