rdf_orig.f90 6.0 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209
  1. module readdata
  2. contains
  3. subroutine readdcd(maxframes,maxatoms,x,y,z,xbox,ybox,zbox,natoms,nframes)
  4. use cudafor
  5. implicit none
  6. integer i,j
  7. integer maxframes,maxatoms
  8. double precision d(6),xbox,ybox,zbox
  9. real*4, managed, allocatable :: x(:,:)
  10. real*4, managed, allocatable :: y(:,:)
  11. real*4, managed, allocatable :: z(:,:)
  12. real*4 dummyr
  13. integer*4 nset, natoms, dummyi,nframes,tframes
  14. character*4 dummyc
  15. open(10,file='/home/hpclabs/manish/input/alk.traj.dcd',status='old',form='unformatted')
  16. read(10) dummyc, tframes,(dummyi,i=1,8),dummyr, (dummyi,i=1,9)
  17. read(10) dummyi, dummyr,dummyr
  18. read(10) natoms
  19. print*,"Total number of frames and atoms are",tframes,natoms
  20. allocate ( x(maxframes,natoms) )
  21. allocate ( y(maxframes,natoms) )
  22. allocate ( z(maxframes,natoms) )
  23. do i = 1,nframes
  24. read(10) (d(j),j=1, 6)
  25. read(10) (x(i,j),j=1,natoms)
  26. read(10) (y(i,j),j=1,natoms)
  27. read(10) (z(i,j),j=1,natoms)
  28. end do
  29. xbox=d(1)
  30. ybox=d(3)
  31. zbox=d(6)
  32. print*,"File reading is done: xbox,ybox,zbox",xbox,ybox,zbox
  33. return
  34. end subroutine readdcd
  35. attributes(global) subroutine pair_calculation( x,y,z,g,natoms,nframes,xbox,ybox,zbox,del,cut)
  36. use cudafor
  37. implicit none
  38. real*4 :: x(:,:)
  39. real*4 :: y(:,:)
  40. real*4 :: z(:,:)
  41. double precision,intent(inout) :: g(:)
  42. integer, value :: nframes,natoms,ind
  43. double precision, value :: xbox,ybox,zbox,del,cut
  44. integer i,j,iconf
  45. double precision dx,dy,dz,r,oldvalue
  46. i = (blockIdx%x-1)*blockDim%x+threadIdx%x
  47. j = (blockIdx%y-1)*blockDim%y+threadIdx%y
  48. if ( i == 1 .and. j == 1) then
  49. print *, natoms,nframes,xbox,ybox,zbox,del,cut, x(1,1), y(1,1), z(1,1), g(1)
  50. end if
  51. do iconf=1,nframes
  52. if(i<=natoms .and. j<=natoms) then
  53. dx=x(iconf,i)-x(iconf,j)
  54. dy=y(iconf,i)-y(iconf,j)
  55. dz=z(iconf,i)-z(iconf,j)
  56. dx=dx-nint(dx/xbox)*xbox
  57. dy=dy-nint(dy/ybox)*ybox
  58. dz=dz-nint(dz/zbox)*zbox
  59. r=dsqrt(dx**2+dy**2+dz**2)
  60. ind=int(r/del)+1
  61. if(r<cut)then
  62. oldvalue = atomicadd(g(ind),1.0d0)
  63. endif
  64. endif
  65. end do
  66. return
  67. end subroutine pair_calculation
  68. end module readdata
  69. program rdf
  70. use readdata
  71. use cudafor
  72. implicit none
  73. integer n,i,j,iconf,ind,istat
  74. integer natoms,nframes,nbin
  75. integer maxframes,maxatoms
  76. type(dim3) :: blocks, threads
  77. parameter (maxframes=10,maxatoms=60000,nbin=2000)
  78. real*4, managed, allocatable :: x(:,:)
  79. real*4, managed, allocatable :: y(:,:)
  80. real*4, managed, allocatable :: z(:,:)
  81. double precision dx,dy,dz
  82. double precision xbox,ybox,zbox,cut
  83. double precision vol,r,del,s2,s2bond
  84. double precision, managed, allocatable :: g(:)
  85. double precision rho,gr,lngr,lngrbond,pi,const,nideal,rf
  86. double precision rlower,rupper
  87. character atmnm*4
  88. real*4 start,finish
  89. open(23,file='RDF.dat',status='unknown')
  90. open(24,file='Pair_entropy.dat',status='unknown')
  91. nframes=10
  92. call cpu_time(start)
  93. print*,"Going to read coordinates"
  94. call readdcd(maxframes,maxatoms,x,y,z,xbox,ybox,zbox,natoms,nframes)
  95. allocate ( g(nbin) )
  96. do i=1,nbin
  97. g(i) = 0.0d0
  98. end do
  99. pi=dacos(-1.0d0)
  100. vol=xbox*ybox*zbox
  101. rho=dble(natoms)/vol
  102. del=xbox/dble(2.0*nbin)
  103. write(*,*) "bin width is : ",del
  104. cut = dble(xbox * 0.5);
  105. threads = dim3(16,16,1)
  106. blocks = dim3(ceiling(real(natoms)/threads%x), ceiling(real(natoms)/threads%y),1)
  107. print *, "natoms:, Grid Size:", threads, blocks
  108. call pair_calculation<<<blocks,threads>>>(x,y,z,g,natoms,nframes,xbox,ybox,zbox,del,cut)
  109. istat = cudaDeviceSynchronize()
  110. if(istat /= 0) then
  111. print *, "Error"
  112. end if
  113. !do iconf=1,nframes
  114. ! do i=1,natoms
  115. ! do j=1,natoms
  116. ! if ( i == 1 .and. j == 1) then
  117. ! print *, natoms,nframes,xbox,ybox,zbox,del,cut, x(1,1), y(1,1), z(1,1), g(1)
  118. ! end if
  119. ! dx=x(iconf,i)-x(iconf,j)
  120. ! dy=y(iconf,i)-y(iconf,j)
  121. ! dz=z(iconf,i)-z(iconf,j)
  122. ! dx=dx-nint(dx/xbox)*xbox
  123. ! dy=dy-nint(dy/ybox)*ybox
  124. ! dz=dz-nint(dz/zbox)*zbox
  125. ! r=dsqrt(dx**2+dy**2+dz**2)
  126. ! ind=int(r/del)+1
  127. ! !if (ind.le.nbin) then
  128. ! if(r<cut)then
  129. ! g(ind)=g(ind)+1.0d0
  130. ! endif
  131. ! enddo
  132. ! enddo
  133. ! enddo
  134. s2=0.01d0
  135. s2bond=0.01d0
  136. const=(4.0d0/3.0d0)*pi*rho
  137. do i=1,nbin
  138. rlower=dble((i-1)*del)
  139. rupper=rlower+del
  140. nideal=const*(rupper**3-rlower**3)
  141. g(i)=g(i)/(dble(nframes)*dble(natoms)*nideal)
  142. r=dble(i)*del
  143. !print*,i+del,i,rupper,rlower,r
  144. !print*,r
  145. if (r.lt.2.0) then
  146. gr=0.0
  147. else
  148. gr=g(i)
  149. endif
  150. if (gr.lt.1e-5) then
  151. lngr=0.0
  152. else
  153. lngr=dlog(gr)
  154. endif
  155. if (g(i).lt.1e-6) then
  156. lngrbond=0.01
  157. else
  158. lngrbond=dlog(g(i))
  159. endif
  160. s2=s2-2*pi*rho*((gr*lngr)-gr+1)*del*r**2.0
  161. s2bond=s2bond-2*pi*rho*((g(i)*lngrbond)-g(i)+1)*del*r*r
  162. rf=dble(i-.5)*del
  163. write(23,*) rf,g(i)
  164. enddo
  165. write(24,*)"s2 : ",s2
  166. write(24,*)"s2bond : ",s2bond
  167. call cpu_time(finish)
  168. print*,"starting at time",start,"and ending at",finish
  169. stop
  170. deallocate(x,y,z,g)
  171. end