rdf_unified_memory.cu 5.1 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177
  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='../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. do iconf=1,nframes
  49. if(i<=natoms .and. j<=natoms) then
  50. dx=x(iconf,i)-x(iconf,j)
  51. dy=y(iconf,i)-y(iconf,j)
  52. dz=z(iconf,i)-z(iconf,j)
  53. dx=dx-nint(dx/xbox)*xbox
  54. dy=dy-nint(dy/ybox)*ybox
  55. dz=dz-nint(dz/zbox)*zbox
  56. r=dsqrt(dx**2+dy**2+dz**2)
  57. ind=int(r/del)+1
  58. if(r<cut)then
  59. oldvalue = atomicadd(g(ind),1.0d0)
  60. endif
  61. endif
  62. end do
  63. return
  64. end subroutine pair_calculation
  65. end module readdata
  66. program rdf
  67. use readdata
  68. use cudafor
  69. implicit none
  70. integer n,i,j,iconf,ind,istat
  71. integer natoms,nframes,nbin
  72. integer maxframes,maxatoms
  73. type(dim3) :: blocks, threads
  74. parameter (maxframes=10,maxatoms=60000,nbin=2000)
  75. real*4, managed, allocatable :: x(:,:)
  76. real*4, managed, allocatable :: y(:,:)
  77. real*4, managed, allocatable :: z(:,:)
  78. double precision dx,dy,dz
  79. double precision xbox,ybox,zbox,cut
  80. double precision vol,r,del,s2,s2bond
  81. double precision, managed, allocatable :: g(:)
  82. double precision rho,gr,lngr,lngrbond,pi,const,nideal,rf
  83. double precision rlower,rupper
  84. character atmnm*4
  85. real*4 start,finish
  86. open(23,file='RDF.dat',status='unknown')
  87. open(24,file='Pair_entropy.dat',status='unknown')
  88. nframes=10
  89. call cpu_time(start)
  90. print*,"Going to read coordinates"
  91. call readdcd(maxframes,maxatoms,x,y,z,xbox,ybox,zbox,natoms,nframes)
  92. allocate ( g(nbin) )
  93. do i=1,nbin
  94. g(i) = 0.0d0
  95. end do
  96. pi=dacos(-1.0d0)
  97. vol=xbox*ybox*zbox
  98. rho=dble(natoms)/vol
  99. del=xbox/dble(2.0*nbin)
  100. write(*,*) "bin width is : ",del
  101. cut = dble(xbox * 0.5);
  102. threads = dim3(16,16,1)
  103. blocks = dim3(ceiling(real(natoms)/threads%x), ceiling(real(natoms)/threads%y),1)
  104. print *, "natoms:, Grid Size:", threads, blocks
  105. call pair_calculation<<<blocks,threads>>>(x,y,z,g,natoms,nframes,xbox,ybox,zbox,del,cut)
  106. istat = cudaDeviceSynchronize()
  107. if(istat /= 0) then
  108. print *, "Error"
  109. end if
  110. s2=0.01d0
  111. s2bond=0.01d0
  112. const=(4.0d0/3.0d0)*pi*rho
  113. do i=1,nbin
  114. rlower=dble((i-1)*del)
  115. rupper=rlower+del
  116. nideal=const*(rupper**3-rlower**3)
  117. g(i)=g(i)/(dble(nframes)*dble(natoms)*nideal)
  118. r=dble(i)*del
  119. if (r.lt.2.0) then
  120. gr=0.0
  121. else
  122. gr=g(i)
  123. endif
  124. if (gr.lt.1e-5) then
  125. lngr=0.0
  126. else
  127. lngr=dlog(gr)
  128. endif
  129. if (g(i).lt.1e-6) then
  130. lngrbond=0.01
  131. else
  132. lngrbond=dlog(g(i))
  133. endif
  134. s2=s2-2*pi*rho*((gr*lngr)-gr+1)*del*r**2.0
  135. s2bond=s2bond-2*pi*rho*((g(i)*lngrbond)-g(i)+1)*del*r*r
  136. rf=dble(i-.5)*del
  137. write(23,*) rf,g(i)
  138. enddo
  139. write(24,*)"s2 : ",s2
  140. write(24,*)"s2bond : ",s2bond
  141. call cpu_time(finish)
  142. print*,"starting at time",start,"and ending at",finish
  143. stop
  144. deallocate(x,y,z,g)
  145. end