-
Notifications
You must be signed in to change notification settings - Fork 1
/
Copy pathoptions_main.f90
602 lines (433 loc) · 22.8 KB
/
options_main.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
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
module options_main
use read_write
use pair_dist
use mod_angles
use mod_check_min
contains
subroutine get_mat_rcut_neighbor(inunit, outunit, natoms, N_species, bins, Rmax, N_file, mat_neighbor, mat_rcut, mat_pdf)
implicit none
integer, intent(in) :: inunit, outunit, natoms, N_species, bins
real*8, intent(in) :: Rmax
real*8, intent(inout) :: mat_rcut(:,:)
integer, intent(inout) ::mat_neighbor(:,:)
integer, intent(out) :: N_file
real*8, intent(out) :: mat_pdf(:,:,:)
real*8, allocatable :: dist_matrix(:,:)
integer, allocatable :: dist_atoms(:,:), N(:)
integer :: numIons(N_species), atomtype(natoms)
real*8 :: pdf(bins), &
cell(3,3), coor(natoms,3)
real*8 :: V, aux, V_mean
integer :: io, type1, type2, min_pdf
character(4) :: caux1, caux2
write(outunit,*) "1. READ TRAJECTORY AND COMPUTE RADIAL DISTRIBUTION FUNCTIONS"
N_file = 0
V_mean = 0
mat_pdf = 0.0d0
do
call read_trj(inunit,cell,coor,numIons,atomtype,io)
if (io < 0) exit
!if(N_file == 11) exit
N_file = N_file + 1
write(outunit,*) "READING FILE ", N_file
call makeMatrices(cell,coor,numIons,atomType,Rmax,N,V,dist_matrix,dist_atoms)
V_mean = V_mean + V
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)*V_mean/N_file
enddo
enddo
rewind(unit=inunit)
!print*, N_file
!print*, "kaka"
! ------------------------------------------------------------
! 2. FIND THE MINIMUM OF THE RDF FOR EACH PAIR OF SPECIES
! ------------------------------------------------------------
write(unit=outunit,fmt=*) "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"//trim(adjustl(caux1))//trim(adjustl(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
aux = integrate(mat_pdf(type1,type2,:),min_pdf,numions(type2)/V,rmax/real(bins))
mat_neighbor(type1,type2) = nint(aux)
write(outunit,"(2i5,2f10.6)") type1, type2, mat_rcut(type1,type2), aux
call N_neighbor_distance("num_neigh_dist"//trim(adjustl(caux1))//trim(adjustl(caux2))//".dat", &
mat_pdf(type1,type2,:), numions(type2)/V, rmax)
enddo
enddo
end subroutine get_mat_rcut_neighbor
subroutine get_N_file(inunit, outunit, natoms, N_file)
implicit none
integer, intent(in) :: inunit, outunit, natoms
integer, intent(out) :: N_file
integer :: io, i
N_file = 0
do
read(inunit, fmt=*, iostat=io)
if (io < 0) exit
do i = 1, 4
read(inunit,fmt=*)
enddo
do i = 1, 3
read(inunit,fmt=*)
enddo
read(inunit,fmt=*)
do i = 1, natoms
read(inunit,fmt=*)
enddo
N_file = N_file + 1
enddo
rewind(unit=inunit)
end subroutine
subroutine get_neighbor_tags(inunit, outunit, natoms, N_species, bins, ext, rmax, mat_neighbor, &
neighbor_order_list, mat_pdf, tot_pdf, contribution_pdf, cont_pdf, &
mat_adf, tot_adf, contribution_adf, cont_adf, &
numIons, atomtype )
implicit none
integer, intent(in) :: inunit, outunit, natoms, N_species, bins, ext, mat_neighbor(:,:)
real*8, intent(in) :: rmax
real*8, intent(out), allocatable :: mat_pdf(:,:,:,:), tot_pdf(:,:,:), contribution_pdf(:,:,:,:), &
mat_adf(:,:,:,:), tot_adf(:,:,:), contribution_adf(:,:,:,:)
integer, intent(out), allocatable :: neighbor_order_list(:,:,:), cont_pdf(:,:,:), cont_adf(:,:,:)
integer, intent(out) :: numIons(N_species), atomtype(natoms)
real*8, allocatable :: dist_matrix(:,:)
integer, allocatable :: dist_atoms(:,:), N(:), neighbor_list(:,:,:,:)
integer :: N_neighbor(natoms, N_species), io
real*8 :: pdf(bins), cell(3,3), coor(natoms,3)
real*8 :: V
integer :: max_neigh, max_pair, i, j, type1, type2
write(outunit,*) "3. READ FIRST CONFIGURATION AND GET THE NEIGHBOR TAGS"
call read_trj(inunit,cell,coor,numIons,atomtype,io)
call makeMatrices(cell,coor,numIons,atomType,Rmax,N,V,dist_matrix,dist_atoms)
call get_neighbor_list2( ext, dist_matrix, dist_atoms, natoms, N_species, atomtype, &
mat_neighbor,N_neighbor, neighbor_list )
allocate(neighbor_order_list(natoms,N_species,maxval(mat_neighbor)+ext))
do i = 1, natoms
do type1 = 1, N_species
do j = 1, 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_neigh = maxval(mat_neighbor)!+ext
max_pair = max_neigh*(max_neigh-1)/2! + ext
max_neigh = max_neigh + ext
allocate( mat_adf(natoms,N_species,max_pair,bins), mat_pdf(natoms,N_species,max_neigh,bins), &
tot_pdf(N_species, N_species, bins), contribution_pdf(N_species,N_species,max_neigh,bins), &
cont_pdf(natoms,N_species,max_neigh), cont_adf(natoms,N_species,max_pair), &
tot_adf(N_species, N_species, bins), contribution_adf(N_species,N_species,max_pair,bins) )
mat_adf = 0.0d0
cont_pdf = 0
cont_adf = 0
contribution_pdf = 0.0d0
contribution_adf = 0.0d0
rewind(unit = inunit)
end subroutine get_neighbor_tags
subroutine get_distributions_dist_angles(inunit, outunit, N_file, natoms, N_species, bins, ext, rmax, mat_neighbor, &
neighbor_order_list, mat_pdf, tot_pdf, contribution_pdf, cont_pdf, &
mat_adf, tot_adf, contribution_adf, cont_adf, V_mean )
integer, intent(in) :: inunit, outunit, natoms, N_species, bins, ext, N_file, mat_neighbor(:,:)
real*8, intent(in) :: rmax
real*8, intent(inout) :: mat_pdf(:,:,:,:), tot_pdf(:,:,:), contribution_pdf(:,:,:,:), &
mat_adf(:,:,:,:), tot_adf(:,:,:), contribution_adf(:,:,:,:)
integer, intent(inout) :: neighbor_order_list(:,:,:), cont_pdf(:,:,:), cont_adf(:,:,:)
real*8, intent(out) :: V_mean
real*8, allocatable :: dist_matrix(:,:)
integer, allocatable :: dist_atoms(:,:), N(:), neighbor_list(:,:,:,:)
integer :: numIons(N_species), atomtype(natoms), N_neighbor(natoms, N_species), io
real*8 :: cell(3,3), coor(natoms,3)
real*8 :: V
integer :: ifile, iat, iesp, N_neigh, N_pair
write(outunit,*) "4. READ TRAJECTORY AND COMPUTE THE DISTRIBUTION OF EACH ANGLE AND EACH BOND"
V_mean = 0.0d0
do ifile = 1, N_file
call read_trj(inunit,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( ext, dist_matrix, dist_atoms, natoms, N_species, atomtype, &
mat_neighbor, N_neighbor, neighbor_list )
V_mean = V_mean + V/N_file
call update_angle_distr ( ext, neighbor_order_list, atomtype, mat_neighbor, dist_matrix, &
neighbor_list, mat_adf, cont_adf )
call update_dist_distr ( ext, V, neighbor_order_list, atomtype, numIons, mat_neighbor, &
dist_matrix, neighbor_list, rmax, mat_pdf, cont_pdf )
deallocate( dist_matrix, dist_atoms, neighbor_list )
enddo
close(unit = inunit)
! Normilize and apply smearing
do iat = 1, natoms
! Loop in species
do iesp = 1, N_species
if ( atomtype(iat) == iesp ) cycle
N_neigh = mat_neighbor(atomtype(iat),iesp)
N_pair = N_neigh*(N_neigh-1)/2
do n1 = 1, N_neigh + ext
call apply_smearing_dist(mat_pdf(iat,iesp,n1,:), rmax, cont_pdf(iat,iesp,n1))
!call normalize_gdr_dist(mat_pdf(iat,iesp,n1,:), rmax, V_mean,numions(iesp), numIons(atomtype(iat)))
enddo
do n1 = 1, N_pair !+ ext
call apply_smearing_ang(mat_adf(iat,iesp,n1,:))
enddo
enddo
enddo
end subroutine get_distributions_dist_angles
subroutine get_deviation_each_dist_angle(outunit, natoms, N_species, atomtype, numions, rmax, V_mean, ext, mat_neighbor, &
mat_pdf, mat_adf, sigma_pdf, sigma_adf)
implicit none
integer, intent(in) :: outunit, natoms, N_species, atomtype(:), numions(:), ext, mat_neighbor(:,:)
real*8, intent(in) :: rmax, V_mean
real*8, intent(inout) :: mat_adf(:,:,:,:), mat_pdf(:,:,:,:)
real*8, intent(out), allocatable :: sigma_pdf(:,:,:), sigma_adf(:,:,:)
real*8, allocatable :: mean_pdf(:,:,:), mean_adf(:,:,:)
integer :: iat, iesp, max_pair, max_neigh, N_neigh, N_pair, jesp
character(4) :: caux1, caux2
write(outunit,*) "5. COMPUTE THE STANDAR DEVIATION OF EACH RADIAL AND ANGLE DISTRIBUTION"
max_neigh = size(mat_pdf,3)
max_pair = size(mat_adf,3)
allocate( mean_adf(natoms,N_species,max_pair), sigma_adf(natoms,N_species,max_pair), &
mean_pdf(natoms,N_species,max_neigh), sigma_pdf(natoms,N_species,max_neigh) )
open(unit = 1, action = "write", status = "replace", file = "deviations.log")
write(1,"(2i5)") natoms, N_species
write(1,"(1000i5)") atomType(:)
do iesp = 1, N_species
write(1,"(100i5)") mat_neighbor(iesp,:)
enddo
write(1,*)
do iat = 1, natoms
do iesp = 1, N_species
N_neigh = mat_neighbor(atomtype(iat),iesp)
N_pair = N_neigh*(N_neigh-1)/2
if ( atomtype(iat) == iesp ) cycle
! N_pair + ext
call get_mean_sigma_angle( N_pair, mat_adf(iat,iesp,:,:), mean_adf(iat,iesp,:), sigma_adf(iat,iesp,:))
call get_mean_sigma_dist( N_neigh, mat_pdf(iat,iesp,:,:), mean_pdf(iat,iesp,:), sigma_pdf(iat,iesp,:), &
rmax, numions(iesp)/V_mean, .true. )
!if (atomType(iat)==2 .and. iesp==1) cycle
write(1,"(2i5,100f10.3)") iat,iesp, mean_adf(iat,iesp,:N_pair)
write(1,"(2i5,100f10.3)") iat,iesp, sigma_adf(iat,iesp,:N_pair)
write(1,"(2i5,100f10.3)") iat,iesp, mean_pdf(iat,iesp,:N_neigh)
write(1,"(2i5,100f10.3)") iat,iesp, sigma_pdf(iat,iesp,:N_neigh)
write(1,*)
write(caux1,"(i0)") iat
write(caux2,"(i0)") iesp
!call write_angle_distr_full( "adf_at"//trim(adjustl(caux1))//"_esp"//trim(adjustl(caux2))//".dat", &
! N_pair, mat_adf(iat,iesp,:,:) )
enddo
enddo
close(unit = 1)
end subroutine get_deviation_each_dist_angle
subroutine get_deviation_contribution(outunit, natoms, N_species, numions, rmax, V_mean, ext, mat_neighbor, &
contribution_pdf, contribution_adf)
implicit none
integer, intent(in) :: outunit, N_species, natoms, numions(:), ext, mat_neighbor(:,:)
real*8, intent(in) :: rmax, V_mean
real*8, intent(inout) :: contribution_pdf(:,:,:,:), contribution_adf(:,:,:,:)
real*8, allocatable :: mean_adf(:,:,:), sigma_adf(:,:,:), mean_pdf(:,:,:), sigma_pdf(:,:,:)
integer :: type1, type2, max_neigh, max_pair, N_pair, N_neigh, n1
character(17) :: datafile
character(4) :: caux1, caux2
max_neigh = size(contribution_pdf,3)
max_pair = size(contribution_adf,3)
allocate( mean_adf(natoms,N_species,max_pair), sigma_adf(natoms,N_species,max_pair), &
mean_pdf(natoms,N_species,max_neigh), sigma_pdf(natoms,N_species,max_neigh) )
do type1 = 1, N_species
do type2 = 1, N_species
if (type1 == type2) cycle
write(caux1,"(i1)") type1
write(caux2,"(i1)") type2
N_neigh = mat_neighbor(type1,type2)
N_pair = N_neigh*(N_neigh-1)/2
call get_mean_sigma_angle( N_pair, contribution_adf(type1,type2,:,:), mean_adf(type1,type2,:), &
sigma_adf(type1,type2,:))
call get_mean_sigma_dist( N_neigh, contribution_pdf(type1,type2,:,:), mean_pdf(type1,type2,:), &
sigma_pdf(type1,type2,:), rmax, numions(type2)/V_mean, .true. )
datafile = "sigma_angle"//trim(adjustl(caux1))//trim(adjustl(caux2))//".dat"
call write_neigh_sigma( datafile, ext, N_pair, mean_adf(type1,type2,:N_pair), &
sigma_adf(type1,type2,:N_pair), .false. )
datafile = "sigma_dista"//trim(adjustl(caux1))//trim(adjustl(caux2))//".dat"
! Rescale sigma --> sigma/mean
! do n1 = 1, N_neigh
! sigma_pdf(type1,type2,n1) = sigma_pdf(type1,type2,n1)/mean_pdf(type1,type2,n1)
! enddo
call write_neigh_sigma( datafile, ext, N_neigh-1, mean_pdf(type1,type2,:N_neigh), &
sigma_pdf(type1,type2,:N_pair), .true. )
enddo
enddo
end subroutine
subroutine dist_distr2pdf(ext, natoms, N_species, mat_neighbor, rmax, V_mean, numions, atomtype, contribution_pdf, tot_pdf)
implicit none
integer, intent(in) :: ext, natoms, N_species, mat_neighbor(:,:), numions(:), atomtype(:)
real*8, intent(in) :: rmax, V_mean
real*8, intent(inout) :: tot_pdf(:,:,:), contribution_pdf(:,:,:,:)
integer :: type1, type2, n1, N_neigh
do type1 = 1, N_species
do type2 = 1, N_species
if (type1 == type2) cycle
N_neigh = mat_neighbor(type1,type2)
do n1 = 1, N_neigh + ext
call normalize_gdr_dist(contribution_pdf(type1,type2,n1,:), rmax, V_mean,numions(type1), numIons(type2))
enddo
call normalize_gdr_dist(tot_pdf(type1,type2,:), rmax, V_mean,numions(type1), numIons(type2))
enddo
enddo
end subroutine dist_distr2pdf
subroutine write_plot_contribution_total (outunit, N_species, mat_neighbor, ext, rmax, V_mean, numions, plot_results, &
contribution_pdf, tot_pdf, contribution_adf, tot_adf)
implicit none
integer, intent(in) :: outunit, N_species, ext, mat_neighbor(:,:), numions(:)
real*8, intent(in) :: rmax, V_mean
real*8, intent(inout) :: tot_pdf(:,:,:), contribution_pdf(:,:,:,:), &
tot_adf(:,:,:), contribution_adf(:,:,:,:)
logical, intent(in) :: plot_results
integer :: type1, type2
character(4) :: caux1, caux2, caux3
character(19) :: datafile
character(17) :: datafile2
real*8, allocatable :: mean(:,:,:), sigma(:,:,:)
integer :: N_pair, max_pair
max_pair = size(contribution_adf,3)
allocate( mean(N_species,N_species,max_pair), sigma(N_species,N_species,max_pair))
open(unit = 11, action = "write", status = "replace", file = "instructions.gnuplot")
do type1 = 1, N_species
do type2 = 1, N_species
if (type1 == type2) cycle
write(caux1,"(i1)") type1
write(caux2,"(i1)") type2
! DISTANCES: DISTRIBUTIONS
!-------------------------------------------------------------------------------
datafile = "contr_tot_gdr"//trim(adjustl(caux1))//trim(adjustl(caux2))//".dat"
call write_gdr_tot_contr( datafile, &
mat_neighbor(type1,type2)+ext, contribution_pdf(type1,type2,:,:), &
tot_pdf(type1,type2,:), rmax )
write(caux3,"(i0)") mat_neighbor(type1,type2)+ext
write(unit = 11,fmt=*) "gnuplot -e "//'"'//"datafile='"//datafile//"'; outfile='"//&
"fig_pdf"//trim(adjustl(caux1))//trim(adjustl(caux2))//".pdf'; nneigh="//&
trim(adjustl(caux3))//'"'//" plot_contr_pdf.gnuplot"
! Gnuplot
if (plot_results) then
call system( "gnuplot -e "//'"'//"datafile='"//datafile//"'; outfile='"//&
"fig_pdf"//trim(adjustl(caux1))//trim(adjustl(caux2))//".pdf'; nneigh="//&
trim(adjustl(caux3))//'"'//" plot_contr_pdf.gnuplot" )
endif
!-------------------------------------------------------------------------------
! DISTANCES: SIGMA
!-------------------------------------------------------------------------------
N_pair = mat_neighbor(type1,type2)!+ext
datafile2 = "sigma_dista"//trim(adjustl(caux1))//trim(adjustl(caux2))//".dat"
write(unit = 11,fmt=*) "gnuplot -e "//'"'//"datafile='"//datafile2//"'; outfile='"//&
"fig_sigma_distance"//trim(adjustl(caux1))//trim(adjustl(caux2))//".pdf'; nneigh="//&
trim(adjustl(caux3))//'"'//" plot_sigma_neigh.gnuplot"
if (plot_results) then
call system( "gnuplot -e "//'"'//"datafile='"//datafile2//"'; outfile='"//&
"fig_sigma_distance"//trim(adjustl(caux1))//trim(adjustl(caux2))//".pdf'; nneigh="//&
trim(adjustl(caux3))//'"'//" plot_sigma_neigh.gnuplot" )
endif
! write(*,"(2i5,100f10.3)") type1, type2, mean(type1,type2,:N_pair)
! write(*,"(2i5,100f10.3)") type1, type2, sigma(type1,type2,:N_pair)
! write(*,*)
!-------------------------------------------------------------------------------
! ANGLES: DISTRIBUTIONS
!-------------------------------------------------------------------------------
datafile = "contr_tot_adf"//trim(adjustl(caux1))//trim(adjustl(caux2))//".dat"
call write_angle_tot_contr( datafile, ext, &
mat_neighbor(type1,type2), contribution_adf(type1,type2,:,:), &
tot_adf(type1,type2,:) )
write(caux3,"(i0)") mat_neighbor(type1,type2)*(mat_neighbor(type1,type2)-1)/2!+ext
write(unit = 11,fmt=*) "gnuplot -e "//'"'//"datafile='"//datafile//"'; outfile='"//&
"fig_adf"//trim(adjustl(caux1))//trim(adjustl(caux2))//".pdf'; nneigh="//&
trim(adjustl(caux3))//'"'//" plot_contr_adf.gnuplot"
if (plot_results) then
call system( "gnuplot -e "//'"'//"datafile='"//datafile//"'; outfile='"//&
"fig_adf"//trim(adjustl(caux1))//trim(adjustl(caux2))//".pdf'; nneigh="//&
trim(adjustl(caux3))//'"'//" plot_contr_adf.gnuplot" )
endif
! ANGLES: SIGMA
!-------------------------------------------------------------------------------
N_pair = mat_neighbor(type1,type2)*(mat_neighbor(type1,type2)-1)/2!+ext
datafile2 = "sigma_angle"//trim(adjustl(caux1))//trim(adjustl(caux2))//".dat"
write(unit = 11,fmt=*) "gnuplot -e "//'"'//"datafile='"//datafile2//"'; outfile='"//&
"fig_sigma_angle"//trim(adjustl(caux1))//trim(adjustl(caux2))//".pdf'; nneigh="//&
trim(adjustl(caux3))//'"'//" plot_sigma_neigh.gnuplot"
if (plot_results) then
call system( "gnuplot -e "//'"'//"datafile='"//datafile2//"'; outfile='"//&
"fig_sigma_angle"//trim(adjustl(caux1))//trim(adjustl(caux2))//".pdf'; nneigh="//&
trim(adjustl(caux3))//'"'//" plot_sigma_neigh.gnuplot" )
endif
! write(*,"(2i5,100f10.3)") type1, type2, mean(type1,type2,:N_pair)
! write(*,"(2i5,100f10.3)") type1, type2, sigma(type1,type2,:N_pair)
! write(*,*)
!-------------------------------------------------------------------------------
enddo
enddo
close (unit = 11)
end subroutine write_plot_contribution_total
subroutine compute_number_constraints(max_sigma_angle)
implicit none
real*8, intent(in) :: max_sigma_angle
integer :: iat, iesp, jesp, natoms, N_species, max_neigh, &
max_pair, N_neigh, N_pair, d1, d2, i
integer, allocatable :: mat_neighbor(:,:), atomtype(:), numions(:), N_const_angle(:,:)
real*8, allocatable :: sigma_angle(:)
open(unit = 1, action = "read", status = "old", file = "deviations.log")
read(1,*) natoms, N_species
allocate(mat_neighbor(N_species, N_species), atomtype(natoms), numions(N_species), &
N_const_angle(N_species, N_species))
read(1,*) atomtype(:)
do iesp = 1, N_species
read(1,*) mat_neighbor(iesp,:)
enddo
read(1,*)
max_neigh = maxval(mat_neighbor)
max_pair = max_neigh*(max_neigh-1)/2
allocate(sigma_angle(max_pair))
numions = 0
N_const_angle = 0
do iat = 1, natoms
numions(atomtype(iat)) = numions(atomtype(iat)) + 1
do iesp = 1, N_species
if ( atomtype(iat) == iesp ) cycle
N_neigh = mat_neighbor(atomtype(iat),iesp)
N_pair = N_neigh*(N_neigh-1)/2
read(1,*)
read(1,*) d1, d2, sigma_angle(:)
read(1,*)
do i = 1, N_pair
if (sigma_angle(i) < max_sigma_angle) then
N_const_angle(atomtype(iat),iesp) = N_const_angle(atomtype(iat),iesp) + 1
else
print*, iat, iesp, sigma_angle(i)
endif
enddo
enddo
enddo
do iesp = 1, N_species
do jesp = 1, N_species
if (iesp == jesp) cycle
N_neigh = mat_neighbor(iesp,jesp)
N_pair = N_neigh*(N_neigh-1)/2
!print*, N_const_angle(iesp,jesp), real(numions(iesp),8), real(N_pair,8)
print*, N_const_angle(iesp,jesp) / real(numions(iesp),8) / real(N_pair,8)
enddo
enddo
!print*, N_const_angle!, numions
close(unit = 1)
end subroutine compute_number_constraints
end module options_main