rdf_offload_collapse.f90 4.7 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166
  1. !/////////////////////////////////////////////////////////////////////////////////////////
  2. !// Author: Manish Agarwal and Gourav Shrivastava , IIT Delhi
  3. !/////////////////////////////////////////////////////////////////////////////////////////
  4. ! Copyright (c) 2021 NVIDIA Corporation. All rights reserved.
  5. module readdata
  6. contains
  7. subroutine readdcd(maxframes,maxatoms,x,y,z,xbox,ybox,zbox,natoms,nframes)
  8. integer i,j
  9. integer maxframes,maxatoms
  10. double precision d(6),xbox,ybox,zbox
  11. real*4, allocatable :: x(:,:)
  12. real*4, allocatable :: y(:,:)
  13. real*4, allocatable :: z(:,:)
  14. real*4 dummyr
  15. integer*4 nset, natoms, dummyi,nframes,tframes
  16. character*4 dummyc
  17. open(10,file='../input/alk.traj.dcd',status='old',form='unformatted')
  18. read(10) dummyc, tframes,(dummyi,i=1,8),dummyr, (dummyi,i=1,9)
  19. read(10) dummyi, dummyr,dummyr
  20. read(10) natoms
  21. print*,"Total number of frames and atoms are",tframes,natoms
  22. allocate ( x(maxframes,natoms) )
  23. allocate ( y(maxframes,natoms) )
  24. allocate ( z(maxframes,natoms) )
  25. do i = 1,nframes
  26. read(10) (d(j),j=1, 6)
  27. read(10) (x(i,j),j=1,natoms)
  28. read(10) (y(i,j),j=1,natoms)
  29. read(10) (z(i,j),j=1,natoms)
  30. end do
  31. xbox=d(1)
  32. ybox=d(3)
  33. zbox=d(6)
  34. print*,"File reading is done: xbox,ybox,zbox",xbox,ybox,zbox
  35. return
  36. end subroutine readdcd
  37. end module readdata
  38. program rdf
  39. use readdata
  40. use nvtx
  41. implicit none
  42. integer n,i,j,iconf,ind
  43. integer natoms,nframes,nbin
  44. integer maxframes,maxatoms
  45. parameter (maxframes=10,maxatoms=60000,nbin=2000)
  46. real*4, allocatable :: x(:,:)
  47. real*4, allocatable :: y(:,:)
  48. real*4, allocatable :: z(:,:)
  49. double precision dx,dy,dz
  50. double precision xbox,ybox,zbox,cut
  51. double precision vol,r,del,s2,s2bond
  52. double precision, allocatable :: g(:)
  53. double precision rho,gr,lngr,lngrbond,pi,const,nideal,rf
  54. double precision rlower,rupper
  55. character atmnm*4
  56. real*4 start,finish
  57. open(23,file='RDF.dat',status='unknown')
  58. open(24,file='Pair_entropy.dat',status='unknown')
  59. nframes=10
  60. call cpu_time(start)
  61. print*,"Going to read coordinates"
  62. call nvtxStartRange("Read File")
  63. call readdcd(maxframes,maxatoms,x,y,z,xbox,ybox,zbox,natoms,nframes)
  64. call nvtxEndRange
  65. allocate ( g(nbin) )
  66. g = 0.0d0
  67. pi=dacos(-1.0d0)
  68. vol=xbox*ybox*zbox
  69. rho=dble(natoms)/vol
  70. del=xbox/dble(2.0*nbin)
  71. write(*,*) "bin width is : ",del
  72. cut = dble(xbox * 0.5);
  73. !pair calculation
  74. !$omp target data map(x(:,:), y (:,:), z (:,:), g (:))
  75. call nvtxStartRange("Pair Calculation")
  76. do iconf=1,nframes
  77. if (mod(iconf,1).eq.0) print*,iconf
  78. !$omp target teams distribute parallel do private(dx,dy,dz,r,ind) collapse(2)
  79. do i=1,natoms
  80. do j=1,natoms
  81. dx=x(iconf,i)-x(iconf,j)
  82. dy=y(iconf,i)-y(iconf,j)
  83. dz=z(iconf,i)-z(iconf,j)
  84. dx=dx-nint(dx/xbox)*xbox
  85. dy=dy-nint(dy/ybox)*ybox
  86. dz=dz-nint(dz/zbox)*zbox
  87. r=dsqrt(dx**2+dy**2+dz**2)
  88. ind=int(r/del)+1
  89. !if (ind.le.nbin) then
  90. if(r<cut)then
  91. !$omp atomic
  92. g(ind)=g(ind)+1.0d0
  93. endif
  94. enddo
  95. enddo
  96. enddo
  97. call nvtxEndRange
  98. !$omp end target data
  99. !entropy calculation
  100. s2=0.01d0
  101. s2bond=0.01d0
  102. const=(4.0d0/3.0d0)*pi*rho
  103. call nvtxStartRange("Entropy Calculation")
  104. do i=1,nbin
  105. rlower=dble((i-1)*del)
  106. rupper=rlower+del
  107. nideal=const*(rupper**3-rlower**3)
  108. g(i)=g(i)/(dble(nframes)*dble(natoms)*nideal)
  109. r=dble(i)*del
  110. if (r.lt.2.0) then
  111. gr=0.0
  112. else
  113. gr=g(i)
  114. endif
  115. if (gr.lt.1e-5) then
  116. lngr=0.0
  117. else
  118. lngr=dlog(gr)
  119. endif
  120. if (g(i).lt.1e-6) then
  121. lngrbond=0.01
  122. else
  123. lngrbond=dlog(g(i))
  124. endif
  125. s2=s2-2*pi*rho*((gr*lngr)-gr+1)*del*r**2.0
  126. s2bond=s2bond-2*pi*rho*((g(i)*lngrbond)-g(i)+1)*del*r*r
  127. rf=dble(i-.5)*del
  128. write(23,*) rf,g(i)
  129. enddo
  130. call nvtxEndRange
  131. write(24,*)"s2 : ",s2
  132. write(24,*)"s2bond : ",s2bond
  133. call cpu_time(finish)
  134. print*,"starting at time",start,"and ending at",finish
  135. stop
  136. deallocate(x,y,z,g)
  137. end