src/ness/ness_stack.f

Fortran project RESTRAX, source module src/ness/ness_stack.f.

Source module last modified on Wed, 13 Jul 2005, 16:20;
HTML image of Fortran source automatically generated by for2html on Mon, 29 May 2006, 15:06.




#--------------------------------------------------------
      logical*4 FUNCTION SPINMATCH(s)
#//// SPIN determines the combination of spin states required:
#//// -3  .... down --> down
#//// -1  .... down --> up
#////  0  .... all
#//// +1  .... up --> down
#//// +3  .... up --> up
#--------------------------------------------------------
      implicit none
      INCLUDE 'ness_common.inc'
      real*4 s
      if(spint.eq.0) then
          SPINMATCH=.true.
      else
          SPINMATCH=(nint(s).eq.nint(spint))
      endif
      end     



#--------------------------------------------------------
      SUBROUTINE EVENT_STACK
# Encapsulates procedures handling event stack:
#        entry    kstack_write:
#        entry    kstack_kf:
#        entry    kstack_n:
#        entry    kstack_allocate:
#        entry    kstack_free:
#        entry    kstack_destroy:
#        entry    getqe:
#        entry    setqe:
#        entry    kstack_phi:
#--------------------------------------------------------
      implicit none
      INCLUDE 'const.inc'
      INCLUDE 'inout.inc'
      INCLUDE 'rescal.inc'
      INCLUDE 'ness_common.inc'

      integer*4 nd,ni,mpage
      parameter(nd=mqom)   ! max. number of records
      parameter(ni=10)   ! number of values in a record
      parameter(mpage=128*1024)   ! size of memory block for events (number of records)

      real*8  k2(3),vq(3),wq(3),vf0,vi0,deh,dev,z,qq
      real*4 GASDEV

      logical*4 map0(3)

      integer*4 indx0,item,i,j,k,ierr,km,n,im,it1,it2,item1,item2
      real*8 qe(4),p,kf,p1,p2,si,sf,vki(3),vkf(3),phi
      logical*4 SPINMATCH,needcopy
      integer*4 kmax(0:mdat),ktop,imax,maxaloc
      real*4 xx
      real, allocatable :: kstore(:,:)  ! auxilliary array for rearranging KSTACK
      real, allocatable :: kstack(:,:)  ! stores events
      save kstack,kmax
      data map0/.true.,.true.,.true./
      data ktop,imax,maxaloc/0,0,0/

# KSTACK stores events for multiple simulations at different configurations:
# first index denotes:
# 1..3  <-- KF(3)
# 4..6  <-- Q(3) in Cooper&Nathans coordinates
#    7  <-- Energy transfer
#    8  <-- PI*PF (probability of KI and KF events)
#    9  <-- Spin transfer (2*SI+SF)
#   10  <-- precession phase difference
 
# second index denotes record index
# KMAX(ITEM) ... number of records allocated for ITEM'th simulation
# KTOP ... total number of allocated records
# IMAX ... number of stored simulations

201   format( 'STACK: ',i7, ' records in ',i3, ' sets, ',g10.4, ' MB')

      ierr=0

#//// Entry to write event


      ENTRY KSTACK_WRITE(indx0,item,vki,vkf,p1,p2,si,sf,phi)
#-----------------------------------------------------------------------
# Entry to write Q,omega, which is calculated from input values of KI,KF,...

      ierr=0
      if(item.gt.imax.or.item.le.0) goto 106
      if(indx0.gt.(kmax(item)-kmax(item-1))) goto 105
        vf0=vkf(1)**2+vkf(2)**2+vkf(3)**2
        vi0=vki(1)**2+vki(2)**2+vki(3)**2
        vq(1)=vkf(1)
        vq(3)=vkf(3)-stp.kf
#* transform dKF (difference from nominal KF) to Lab coord.
        k2(3)=-vq(1)*somega+vq(3)*comega
        k2(2)=vkf(2)   
        k2(1)=+vq(1)*comega+vq(3)*somega
#* get VQ (total Q in lab. coordinates) 
        do i=1,2
            vq(i)=k2(i)-vki(i)
        end do
#* VQ is the difference from nominal Q value 
        vq(3)=k2(3)-vki(3)+stp.ki
      
#* transform Q to C&N coord.
        do i=1,3
            wq(i)=0
            do j=1,3
              wq(i)=wq(i)+mlc(j,i)*vq(j)
            enddo
        enddo
#* add sample mosaicity 
        if(smos.ne.0) then
          qq=sqrt(wq(1)**2+wq(2)**2+wq(3)**2)
          deh=stp.q*smos*GASDEV()   
          dev=stp.q*smos*GASDEV()
#          write(*,*) WQ(2), DEH
#          write(*,*) WQ(3), DEV
#          pause
          wq(2)=wq(2)+deh
          wq(3)=wq(3)+dev
          z=sqrt(wq(1)**2+wq(2)**2+wq(3)**2)
          do i=1,3
            wq(i)=wq(i)*qq/z   ! keep |Q| unchanged
          enddo
        endif    
#* store in KSTACK
        do i=1,3
          xx=vkf(i)
          kstack(i,indx0+kmax(item-1))=xx
          xx=wq(i)
          kstack(i+3,indx0+kmax(item-1))=xx
        enddo

        xx=hsqov2m*(vi0-vf0)-stp.e
        kstack(7,indx0+kmax(item-1))=xx
        xx=p1*p2
        kstack(8,indx0+kmax(item-1))=xx
        xx=2*nint(si)+nint(sf)
        kstack(9,indx0+kmax(item-1))=xx
        xx=phi
        kstack(10,indx0+kmax(item-1))=xx
      return


      ENTRY KSTACK_KF(indx0,item,kf,p)
#--------------------------------------------------------------------------
# Entry to read Kf
      ierr=0
      if(item.gt.imax.or.item.le.0) goto 106
      if(indx0.gt.(kmax(item)-kmax(item-1))) goto 105
         kf=0
         do i=1,3
           kf=kf+kstack(i,indx0+kmax(item-1))**2
         enddo
         kf=sqrt(kf)
         if(SPINMATCH(kstack(9,indx0+kmax(item-1)))) then
            p=kstack(8,indx0+kmax(item-1))
         else
            p=0
         endif
      return


      ENTRY KSTACK_N(indx0,item)
#--------------------------------------------------------------------------
# Entry to return number of allocated events
      indx0=kmax(item)-kmax(item-1)
      return


      ENTRY KSTACK_ALLOCATE(indx0,item)
#--------------------------------------------------------------------------
# Entry for memory allocation 
# INDX0=num. of records to be allocated for ITEM-th set
      ierr=0
      if(item.lt.1.or.indx0.eq.(kmax(item)-kmax(item-1))) then    ! no reallocation necessary
         return
      endif

#* DEFINE NUMBER OF NEWLY ALLOCATED ROWS
      if(item.gt.mdat) then
         write(smes,*)  'Maximum number of data sets is ',mdat
         stop
      endif
      im=imax  ! new top item
      if(item.gt.imax) im=item  ! a new item will be added on top
      if((indx0+ktop).gt.nd) then
        indx0=nd-ktop
        write(smes,*)  'Maximum allowed number of events: ',indx0
      endif
      km=kmax(item)-kmax(item-1)    ! old number of records in ITEM-th item
      n=ktop-km+indx0               ! new total number of records (= new  KTOP)

#* ALLOCATE NEW MEMORY if necessary       
      if (n.gt.maxaloc) then
        needcopy=(ktop.ge.1.and.im.gt.1.and.(item.ne.1.or.imax.gt.1))           
        if (needcopy) then   !* save current KSTACK in temporary array
        allocate (kstore(1:ni,1:ktop),stat=ierr)
          if (ierr.ne.0) then
            write(smes,*)  'Cannot allocate more memory'
            goto 99
          endif
          do i=1,ktop
            do j=1,ni
              kstore(j,i)=kstack(j,i)
            enddo
          enddo
        endif
        if (allocated(kstack)) deallocate(kstack,stat=ierr)
        if (ierr.ne.0) goto 98
      i=int(n/mpage)+1
      maxaloc=i*mpage
        allocate (kstack(1:ni,1:maxaloc),stat=ierr)
        write(smes,*)  'Allocated memory for ',maxaloc, ' records.'
      if (ierr.ne.0) goto 99
        if (needcopy) then   !* restore current KSTACK from temporary array
          do i=1,ktop
            do j=1,ni
              kstack(j,i)=kstore(j,i)
            enddo
          enddo
          if (allocated(kstore)) deallocate(kstore,stat=ierr)
          if (ierr.ne.0) goto 98
      endif
      endif
           
#* RE-ARRANGE VALUES IN KSTACK 
      if (item.lt.imax.and.km.ne.indx0) then
        if (indx0.lt.km) then  ! shift down
         do i=kmax(item+1)+1,kmax(imax)
             do j=1,ni
               kstack(j,i+indx0-km)=kstack(j,i)
             enddo
         enddo  
        else if (indx0.gt.km) then  ! shift up
         do i=kmax(imax),kmax(item+1)+1,-1
             do j=1,ni
                kstack(j,i+indx0-km)=kstack(j,i)
             enddo      
         enddo
        endif
      endif
#* UPDATE KMAX(I),IMAX AND KTOP 
      do k=item,mdat
         kmax(k)=kmax(k)+indx0-km
      enddo
      imax=im
      ktop=n
#* set the new records in ITEM = 0       
      do i=kmax(item-1)+1,kmax(item)  
         do j=1,ni
           kstack(j,i)=0.
         enddo      
      enddo   

      return


      ENTRY KSTACK_FREERANGE(item1,item2)
#--------------------------------------------------------------------------
# deallocate range of ITEM1..ITEM2 data sets
      ierr=0
      it1=item1
      it2=item2
      if (item2.gt.imax) it2=imax
      if (item1.lt.2) it1=1
      
      if (it2.lt.it1.or.it1.gt.imax.or.it2.lt.1) return 
      
      km=kmax(it2)-kmax(it1-1)    ! number of records in deleted items
      ktop=ktop-km                ! new total number of records 

#* RE-ARRANGE VALUES IN KSTACK 
      if (it2.lt.imax) then
        do i=kmax(it2+1)+1,kmax(imax)
          do j=1,ni
            kstack(j,i-km)=kstack(j,i)
          enddo
      enddo  
      endif

#* UPDATE KMAX(I),IMAX
      imax=imax-1-it2+it1
      do k=it1,imax
         kmax(k)=kmax(k+1+it2-it1)-km
      enddo
      do k=imax+1,mdat
         kmax(k)=kmax(k-1)
      enddo
      return



      ENTRY KSTACK_DESTROY
#--------------------------------------------------------------------------
# Entry for memory deallocation
      kmax=0
      ierr=0
      if (allocated(kstack)) deallocate (kstack,stat=ierr)
      if (ierr.ne.0) goto 98
      if (allocated(kstore))  deallocate (kstore,stat=ierr)
      if (ierr.ne.0) goto 98
      ktop=0
      imax=0
      maxaloc=0
      return


      ENTRY GETQE(indx0,item,qe,p)
#--------------------------------------------------------------------------
# Entry for receive QE(4) vectors with appropriate weight
      if(item.gt.imax) goto 106
      if(indx0.gt.(kmax(item)-kmax(item-1))) goto 105
      if (SPINMATCH(kstack(9,indx0+kmax(item-1)))) then
         do i=1,4
           qe(i)=kstack(i+3,indx0+kmax(item-1))
         enddo
         p=kstack(8,indx0+kmax(item-1))         
      else
         p=0.
      endif
      return


      ENTRY SETQE(indx0,item,qe,p)
#-----------------------------------------------------------------
# Entry for setting QE(4) vector with weight
      if(item.gt.imax) goto 106
      if(indx0.gt.(kmax(item)-kmax(item-1))) goto 105
      do i=1,4
           kstack(i+3,indx0+kmax(item-1))=qe(i)
      enddo
      kstack(8,indx0+kmax(item-1))=p
      kstack(9,indx0+kmax(item-1))=0.
      return

      ENTRY KSTACK_PHI(indx0,item,phi,p)
#-----------------------------------------------------------------
# Entry to receive precession phase difference with appropriate weight
      if(item.gt.imax) goto 106
      if(indx0.gt.(kmax(item)-kmax(item-1))) goto 105
      phi=kstack(10,indx0+kmax(item-1))
      p=kstack(8,indx0+kmax(item-1))
      return

#-----------------------------------------------------------------
# Error messages

98    write(smes,*) 'Error: deallocating memory for event storage: ',
     *    ierr
      stop
#      RETURN

99    write(smes,*)  'Error: allocating memory for event storage: ',
     *    ierr
      stop
#      RETURN

102   write(smes,*) 'Attempt to write to a nonallocated arrea !'
      stop
#      RETURN
103   write(smes,*) 'Attempt to read from a nonallocated arrea !'
      stop
104   write(smes,*) 'Max. 10 event sets are allowed !'
      stop
105   write(smes,*) 'Not allocated: INDX0=',indx0, ' ITEM=',item,
     *   ' IMIN=',kmax(item-1), ' IMAX=',kmax(item)
      stop
106   write(smes,*) 'Not allocated ITEM: ',item
      stop


#      RETURN



      end