-
Notifications
You must be signed in to change notification settings - Fork 1
/
main_lammps.f90
222 lines (167 loc) · 7.45 KB
/
main_lammps.f90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
program main
use read_write
use pair_dist
use mod_check_min
use mod_angles
implicit none
integer :: natoms, io, bins, i, j, iat, iesp
integer :: type1, type2
integer, dimension(:), allocatable :: numIons, atomtype, aux_numions, aux_atomtype
real*8, dimension(:), allocatable :: pdf, angle_distr
real*8, dimension(3,3) :: cell
real*8, dimension(:,:), allocatable :: coor
real*8, parameter :: pi = acos(-1.0d0)
real*8 :: rmax, V
integer :: N_species, min_pdf, N_file, ifile, max_pair, N_pair
integer, dimension(:), allocatable :: N
integer, dimension(:,:), allocatable :: dist_atoms
real*8, dimension(:,:), allocatable :: dist_matrix
integer, allocatable :: neighbor_list(:,:,:,:), N_neighbor(:,:), &
neighbor_order_list(:,:,:)
real*8, allocatable :: mat_rcut(:,:), mat_theta(:,:,:), mat_pdf(:,:,:), mat_adf(:,:,:,:), &
mean(:,:,:), sigma(:,:,:), mat_neighbor(:,:)
character(len=2), allocatable :: species(:)
character(len=100), allocatable :: file_list(:)
character(4) :: caux1, caux2
character(13) :: infile
logical :: same
integer :: nn, ext, cont
infile = "qe1.lammpstrj"
rmax = 5.0
bins = 200
ext = 2
allocate(species(2))
species = (/ "Si", "O " /)
N_species = size(species)
allocate(aux_numions(N_species))
allocate(mat_rcut(N_species,N_species), mat_neighbor(N_species, N_species), mat_theta(N_species,N_species,2))
allocate(mat_pdf(N_species,N_species,bins))
allocate(pdf(bins), angle_distr(bins))
mat_pdf = 0.0d0
! 1. READ TRAJECTORY AND COMPUTE RADIAL DISTRIBUTION FUNCTIONS
! ------------------------------------------------------------
print*, "1. READ TRAJECTORY AND COMPUTE RADIAL DISTRIBUTION FUNCTIONS"
call atom_number(infile,natoms)
allocate(coor(natoms,3),atomtype(natoms), numions(N_species),N_neighbor(natoms,N_species))
open(unit=123,file=infile,status='old',action='read')
! call read_trj(123,cell,coor,numIons,atomtype,io)
! print*, numIons(:)
! open(unit=124,file="kaka_POSCAR",status='replace',action='write')
! call write_vasp(124,cell,coor)
N_file = 0
do
call read_trj(123,cell,coor,numIons,atomtype,io)
if (io < 0) exit
if(N_file == 10) exit
N_file = N_file + 1
print*, "READING FILE ", N_file
call makeMatrices(cell,coor,numIons,atomType,Rmax,N,V,dist_matrix,dist_atoms)
do type1 = 1, N_species
do type2 = 1, N_species
if (type1 == type2) cycle
! Compute radial distribution function
call compute_gdr(dist_matrix,dist_atoms,atomType,type1,type2,numIons,V,bins,rmax,pdf)
mat_pdf(type1,type2,:) = mat_pdf(type1,type2,:) + pdf(:)
enddo
enddo
deallocate(dist_matrix, dist_atoms )
enddo
do type1 = 1, N_species
do type2 = 1, N_species
if (type1==type2) cycle
mat_pdf(type1,type2,:) = mat_pdf(type1,type2,:)/real(N_file,8)
enddo
enddo
rewind(unit=123)
!print*, N_file
!print*, "kaka"
! ------------------------------------------------------------
! 2. FIND THE MINIMUM OF THE RDF FOR EACH PAIR OF SPECIES
! ------------------------------------------------------------
print*, "2. FIND THE MINIMUM OF THE RDF FOR EACH PAIR OF SPECIES"
do type1 = 1, N_species
do type2 = 1, N_species
if (type1 == type2) cycle
write(caux1,"(i1)") type1
write(caux2,"(i1)") type2
call write_gdr( "gdr"//caux1//caux2//".dat", mat_pdf(type1,type2,:), rmax )
! Find the first minimum
call check_min(mat_pdf(type1,type2,:), int(bins*0.5/rmax), min_pdf)
mat_rcut(type1,type2) = rmax/bins*min_pdf
mat_neighbor(type1,type2) = integrate(mat_pdf(type1,type2,:),min_pdf,numions(type2)/V,rmax/real(bins))
write(*,"(2i5,2f10.6)") type1, type2, mat_rcut(type1,type2), mat_neighbor(type1,type2)
enddo
enddo
! ------------------------------------------------------------
! 3. READ FIRST CONFIGURATION AND GET THE NEIGHBOR TAGS
! ------------------------------------------------------------
print*, "3. READ FIRST CONFIGURATION AND GET THE NEIGHBOR TAGS"
call read_trj(123,cell,coor,numIons,atomtype,io)
call makeMatrices(cell,coor,numIons,atomType,Rmax,N,V,dist_matrix,dist_atoms)
call get_neighbor_list2( dist_matrix, dist_atoms, natoms, N_species, atomtype, ceiling(mat_neighbor),N_neighbor, neighbor_list )
allocate(neighbor_order_list(natoms,N_species,maxval(ceiling(mat_neighbor))+ext))
do i = 1, natoms
do type1 = 1, N_species
do j = 1, ceiling(mat_neighbor(atomtype(i),type1))+ext
neighbor_order_list(i,type1,j) = neighbor_list(i,type1,j,2)
enddo
enddo
enddo
deallocate( dist_matrix, dist_atoms, neighbor_list )
max_pair = maxval(ceiling(mat_neighbor))+ext
max_pair = max_pair*(max_pair-1)/2
allocate( mat_adf(natoms,N_species,max_pair,bins) )
mat_adf = 0.0d0
! ------------------------------------------------------------
! 4. READ TRAJECTORY AND COMPUTE THE DISTRIBUTION OF EACH ANGLE
! ------------------------------------------------------------
print*, "4. READ TRAJECTORY AND COMPUTE THE DISTRIBUTION OF EACH ANGLE"
rewind(unit = 123)
!cont = 0
do ifile = 1, 100
call read_trj(123,cell,coor,numIons,atomtype,io)
if (io < 0) exit
if(ifile == 10) exit
call makeMatrices(cell,coor,numIons,atomType,Rmax,N,V,dist_matrix,dist_atoms)
call get_neighbor_list2( dist_matrix, dist_atoms, natoms, N_species, atomtype, ceiling(mat_neighbor),N_neighbor, neighbor_list )
! do i = 1, natoms
! do type1 = 1, N_species
! if (type1 == atomtype(i)) cycle
! nn = ceiling( mat_neighbor(atomtype(i),type1) )
! call compare( neighbor_order_list(i,type1,:nn+ext), neighbor_list(i,type1,:nn,2), same )
! if ( .not. same ) then
! cont = cont + 1
! print*, ifile,i, type1
! print*, neighbor_order_list(i,type1,:nn+ext)
! print*, neighbor_list(i,type1,:nn,2)
! print*,
! endif
! enddo
! enddo
call update_angle_distr ( ext, neighbor_order_list, atomtype, ceiling(mat_neighbor), dist_matrix, &
neighbor_list, mat_adf )
deallocate( dist_matrix, dist_atoms, neighbor_list )
enddo
close(unit = 123)
!print*, cont
! ------------------------------------------------------------
! 5. COMPUTE THE STANDAR DEVIATION OF EACH ANGLE DISTRIBUTION
! ------------------------------------------------------------
print*, "5. COMPUTE THE STANDAR DEVIATION OF EACH ANGLE DISTRIBUTION"
allocate( mean(natoms,N_species,max_pair), sigma(natoms,N_species,max_pair) )
do iat = 1, natoms
do iesp = 1, N_species
N_pair = mat_neighbor(atomtype(iat),iesp)*(mat_neighbor(atomtype(iat),iesp)-1)/2
if ( atomtype(iat) == iesp ) cycle
call get_mean_sigma( N_pair, mat_adf(iat,iesp,:,:), mean(iat,iesp,:), sigma(iat,iesp,:) )
!write(*,"(2i5,100f10.3)") iat,iesp, mean(iat,iesp,:N_pair)
!write(*,"(2i5,100f10.3)") iat,iesp, sigma(iat,iesp,:N_pair)
!write(*,*)
write(caux1,"(i4)") iat
write(caux2,"(i4)") iesp
call write_angle_distr_full( "adf_at"//trim(adjustl(caux1))//"_esp"//trim(adjustl(caux1))//".dat", &
N_pair, mat_adf(iat,iesp,:,:) )
enddo
enddo
! ------------------------------------------------------------
end program