-
Notifications
You must be signed in to change notification settings - Fork 6
/
Copy pathsample_mod.f90
936 lines (659 loc) · 22 KB
/
sample_mod.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
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
module sample_mod
use global_mod
use system_mod
use pbc_mod
implicit none
contains
!-----------------------------------------------------------------------
subroutine PotentialEnergy(trap,VTable,R,Pot,F2)
implicit none
logical :: trap
real (kind=8) :: Pot
real (kind=8) :: rij2,rij
real (kind=8) :: fij
real (kind=8) :: Interpolate
integer (kind=4) :: ip,jp
integer (kind=4) :: k
real (kind=8),optional :: F2
real (kind=8),dimension(0:Nmax+1) :: VTable
real (kind=8),dimension(dim,Np) :: R,F
real (kind=8),dimension(dim) :: xij
Pot = 0.d0
do ip=1,Np
do k=1,dim
if (trap) then
F(k,ip) = TrapPot(1,a_ho(k),R(k,ip))
Pot = Pot+TrapPot(0,a_ho(k),R(k,ip))
else
F(k,ip) = 0.d0
end if
end do
end do
do ip=1,Np-1
do jp=ip+1,Np
rij2 = 0.d0
do k=1,dim
xij(k) = R(k,ip)-R(k,jp)
end do
if (trap) then
rij2 = 0.d0
do k=1,dim
rij2 = rij2+xij(k)*xij(k)
end do
else
call MinimumImage(xij,rij2)
end if
if (trap) then
rij = sqrt(rij2)
if (v_table) then
Pot = Pot+Interpolate(0,Nmax,dr,VTable,rij)
else
Pot = Pot+Potential(xij,rij)
end if
if (present(F2)) then
if (v_table) then
do k=1,dim
fij = Interpolate(1,Nmax,dr,VTable,rij)*xij(k)/rij
F(k,ip) = F(k,ip)+fij
F(k,jp) = F(k,jp)-fij
end do
else
do k=1,dim
fij = Force(k,xij,rij)
F(k,ip) = F(k,ip)+fij
F(k,jp) = F(k,jp)-fij
end do
end if
end if
else
if (rij2<=rcut2) then
rij = sqrt(rij2)
if (v_table) then
Pot = Pot+Interpolate(0,Nmax,dr,VTable,rij)
else
Pot = Pot+Potential(xij,rij)
end if
if (present(F2)) then
if (v_table) then
do k=1,dim
fij = Interpolate(1,Nmax,dr,VTable,rij)*xij(k)/rij
F(k,ip) = F(k,ip)+fij
F(k,jp) = F(k,jp)-fij
end do
else
do k=1,dim
fij = Force(k,xij,rij)
F(k,ip) = F(k,ip)+fij
F(k,jp) = F(k,jp)-fij
end do
end if
end if
end if
end if
end do
end do
if (present(F2)) then
F2 = 0.d0
do ip=1,Np
do k=1,dim
F2 = F2+F(k,ip)*F(k,ip)
end do
end do
end if
return
end subroutine PotentialEnergy
!-----------------------------------------------------------------------
subroutine LocalEnergy(trap,LogWF,VTable,R,E,Kin,Pot)
implicit none
logical :: trap
real (kind=8) :: E,Kin,Pot
real (kind=8) :: rij,rij2
real (kind=8) :: fij
real (kind=8) :: dudr,d2udr2
real (kind=8) :: LapLogPsi
real (kind=8) :: Interpolate
integer (kind=4) :: i,j
integer (kind=4) :: k
real (kind=8),dimension (0:Nmax+1) :: LogWF,VTable
real (kind=8),dimension (dim,Np) :: R,F
real (kind=8),dimension (dim) :: xij
Kin = 0.d0
Pot = 0.d0
LapLogPsi = 0.d0
do i=1,Np
do k=1,dim
if (trap) then
F(k,i) = TrapPsi(1,a_ho(k),R(k,i))
Pot = Pot+TrapPot(0,a_ho(k),R(k,i))
LapLogPsi = LapLogPsi+TrapPsi(2,a_ho(k),R(k,i))
else
F(k,i) = 0.d0
end if
end do
end do
LapLogPsi = 0.5d0*LapLogPsi
do i=1,Np-1
do j=i+1,Np
!Calculation of the distance between each pair of dipoles
!choosing only the nearest image of each particle
do k=1,dim
xij(k) = R(k,i)-R(k,j)
end do
if (trap) then
rij2 = 0.d0
do k=1,dim
rij2 = rij2+xij(k)*xij(k)
end do
else
call MinimumImage(xij,rij2)
end if
!Calculation of the local kinetic energy
!We use in this routine the following prescription:
! Kin = -Laplacian(Psi)/(2*Psi)
! T = -Laplacian(log(Psi))/4
! F = Grad(Psi)/(sqrt(2)*Psi)
! Kin = 2*T-F*F
!We use linear interpolation between the nearest points of the
!tabulated wave function, as in the previous subroutine we use:
!
! x_{i-1} < z < x_i
!
! psi(z) = psi(x_i)*(z-x_{i-1})/dr+psi(x_{i-1})*(x_i-z)/dr
!
!The evaluation of derivatives is also numeric by interpolating
!in the next and previous intervals.
if (trap) then
rij = sqrt(rij2)
if (wf_table) then
dudr = Interpolate(1,Nmax,dr,LogWF,rij)
d2udr2 = Interpolate(2,Nmax,dr,LogWF,rij)
else
dudr = LogPsi(1,Rm,rij)
d2udr2 = LogPsi(2,Rm,rij)
end if
LapLogPsi = LapLogPsi+((real(dim)-1)*dudr/rij+d2udr2)
do k=1,dim
fij = dudr*xij(k)/rij
F(k,i) = F(k,i)+fij
F(k,j) = F(k,j)-fij
end do
!Calculation of the local potential energy
if (v_table) then
Pot = Pot+Interpolate(0,Nmax,dr,VTable,rij)
else
Pot = Pot+Potential(xij,rij)
end if
else
if (rij2<=rcut2) then
rij = sqrt(rij2)
if (wf_table) then
dudr = Interpolate(1,Nmax,dr,LogWF,rij)
d2udr2 = Interpolate(2,Nmax,dr,LogWF,rij)
else
dudr = LogPsi(1,Rm,rij)
d2udr2 = LogPsi(2,Rm,rij)
end if
LapLogPsi = LapLogPsi+((real(dim)-1)*dudr/rij+d2udr2)
do k=1,dim
fij = dudr*xij(k)/rij
F(k,i) = F(k,i)+fij
F(k,j) = F(k,j)-fij
end do
!Calculation of the local potential energy
if (v_table) then
Pot = Pot+Interpolate(0,Nmax,dr,VTable,rij)
else
Pot = Pot+Potential(xij,rij)
end if
end if
end if
end do
end do
!Kinetic energy and drift force calculation
Kin = 2.d0*LapLogPsi
do i=1,Np
do k=1,dim
Kin = Kin+F(k,i)*F(k,i)
end do
end do
!Computing energy
Kin = -0.5d0*Kin
E = Kin+Pot
return
end subroutine LocalEnergy
!-----------------------------------------------------------------------
subroutine ThermEnergy(trap,VTable,Path,dt,E,Ec,Ep)
implicit none
logical :: trap
real (kind=8) :: E,Ec,Ep
real (kind=8) :: dt
real (kind=8) :: Pot,F2
real (kind=8) :: rij2
integer (kind=4) :: ib
integer (kind=4) :: ip
integer (kind=4) :: k
real (kind=8),dimension(0:Nmax+1) :: VTable
real (kind=8),dimension(dim,Np,0:2*Nb) :: Path
real (kind=8),dimension(dim) :: xij
Ep = 0.d0
Ec = 0.d0
E = 0.d0
do ib=0,2*Nb-1
if (mod(ib,2)==0) then
call PotentialEnergy(trap,VTable,Path(:,:,ib),Pot)
F2 = 0.d0
else
call PotentialEnergy(trap,VTable,Path(:,:,ib),Pot,F2)
end if
if (ib==Nb) then
Ep = Pot
end if
E = E+GreenFunction(1,ib,dt,Pot,F2)
do ip=1,Np
rij2 = 0.d0
do k=1,dim
xij(k) = Path(k,ip,ib)-Path(k,ip,ib+1)
end do
if (trap) then
rij2 = 0.d0
do k=1,dim
rij2 = rij2+xij(k)*xij(k)
end do
E = E-0.5d0*rij2/(dt*dt)
else
call MinimumImage(xij,rij2)
if (rij2<=rcut2) E = E-0.5d0*rij2/(dt*dt)
end if
end do
end do
E = 0.5d0*(E/real(Nb)+real(dim*Np)/dt)
Ec = E-Ep
return
end subroutine ThermEnergy
!-----------------------------------------------------------------------
subroutine PairCorrelation(R,gr)
implicit none
real (kind=8) :: rij2,rij
integer (kind=4) :: ip,jp
integer (kind=4) :: k
integer (kind=4) :: ibin
real (kind=8),dimension (dim) :: xij
real (kind=8),dimension (dim,Np) :: R
real (kind=8),dimension (Nbin) :: gr
do ip=1,Np-1
do jp=ip+1,Np
rij2 = 0.d0
do k=1,dim
xij(k) = R(k,ip)-R(k,jp)
end do
call MinimumImage(xij,rij2)
if (rij2<=rcut2) then
rij = sqrt(rij2)
ibin = int(rij/rbin)+1
gr(ibin) = gr(ibin)+2.d0
end if
end do
end do
return
end subroutine PairCorrelation
!-----------------------------------------------------------------------
subroutine StructureFactor(Nk,R,Sk)
implicit none
real (kind=8) :: qr
real (kind=8) :: SumCos,SumSin
integer (kind=4) :: iq,Nk
integer (kind=4) :: k
integer (kind=4) :: ip
real (kind=8),dimension(dim) :: x
real (kind=8),dimension(dim,Np) :: R
real (kind=8),dimension(dim,Nk) :: Sk
SumCos = 0.d0
SumSin = 0.d0
do iq=1,Nk
do k=1,dim
SumCos = 0.d0
SumSin = 0.d0
do ip=1,Np
x(k) = R(k,ip)
qr = real(iq)*qbin(k)*x(k)
SumCos = SumCos+cos(qr)
SumSin = SumSin+sin(qr)
end do
Sk(k,iq) = Sk(k,iq)+(SumCos*SumCos+SumSin*SumSin)
end do
end do
return
end subroutine StructureFactor
!-----------------------------------------------------------------------
subroutine OBDM(xend,nrho)
implicit none
real (kind=8) :: rij2,rij
real (kind=8) :: costheta,sintheta
real (kind=8) :: cos2mtheta,sin2mtheta
complex (kind=8) :: exptheta,exp2theta,exp2mtheta
integer (kind=4) :: ibin,k,m
real (kind=8), dimension (dim,2) :: xend
real (kind=8), dimension (dim) :: xij
real (kind=8), dimension (0:Npw,Nbin) :: nrho
do k=1,dim
xij(k) = xend(k,1)-xend(k,2)
end do
call MinimumImage(xij,rij2)
if (rij2<=rcut2) then
rij = sqrt(rij2)
ibin = int(rij/rbin)+1
sintheta = xij(2)/rij
costheta = xij(1)/rij
exptheta = cmplx(costheta,sintheta,8)
exp2theta = exptheta*exptheta
exp2mtheta = 1.d0
do m=0,Npw
cos2mtheta = real(exp2mtheta)
sin2mtheta = aimag(exp2mtheta)
nrho(m,ibin) = nrho(m,ibin)+cos2mtheta
exp2mtheta = exp2mtheta*exp2theta
end do
end if
return
end subroutine OBDM
!-----------------------------------------------------------------------
subroutine PermutationSampling(new_perm_cycle,end_perm_cycle,iperm,&
&Particles_in_perm_cycle,Perm_histogram,&
&isopen,iw,ik,swap_accepted)
implicit none
logical :: new_perm_cycle
logical :: end_perm_cycle
logical :: isopen
logical :: particle_already_in_cycle
logical, optional :: swap_accepted
integer (kind=4) :: ip
integer (kind=4) :: iw,iperm
integer (kind=4), optional :: ik
integer (kind=4), dimension (Np) :: Particles_in_perm_cycle
integer (kind=4), dimension (Np) :: Perm_histogram
if (new_perm_cycle) then
Particles_in_perm_cycle(:) = 0
Particles_in_perm_cycle(1) = iw
iperm = 1
new_perm_cycle = .false.
end if
if (present(swap_accepted)) then
if (swap_accepted) then
! Check if particle is already in the permutation cycle
! if it is present then end the permutation cycle
do ip=1,Np
if (Particles_in_perm_cycle(ip)==ik) then
!end_perm_cycle = .true.
particle_already_in_cycle = .true.
exit
else
!end_perm_cycle = .false.
particle_already_in_cycle = .false.
end if
end do
if (end_perm_cycle .eqv. .false.) then
if (.not. particle_already_in_cycle) then
iperm = iperm+1
Particles_in_perm_cycle(iperm) =ik
end if
end if
end if
end if
if (end_perm_cycle) then
Perm_histogram(iperm) = Perm_histogram(iperm)+1
if (isopen) then
Particles_in_perm_cycle(:) = 0
Particles_in_perm_cycle(1) = iw
iperm = 1
end if
end_perm_cycle = .false.
end if
return
end subroutine PermutationSampling
!-----------------------------------------------------------------------
!!$ subroutine DensityProfile(R,dens)
!!$
!!$ implicit none
!!$
!!$ real (kind=8) :: x,y
!!$ integer (kind=4) :: ip,ibin,jbin
!!$
!!$ real (kind=8), dimension (dim,Np) :: R
!!$ real (kind=8), dimension (Nbin,Nbin) :: dens
!!$
!!$ do ip=1,Np
!!$
!!$ x = R(1,ip)
!!$ y = R(2,ip)
!!$
!!$ ibin = int((x+0.5d0*rcut)/rbin)
!!$ jbin = int((y+0.5d0*rcut)/rbin)
!!$
!!$ if (ibin>0) then
!!$ if (ibin<=Nbin) then
!!$ if (jbin>0) then
!!$ if (jbin<=Nbin) then
!!$ dens(ibin,jbin) = dens(ibin,jbin)+1.d0
!!$ end if
!!$ end if
!!$ end if
!!$ end if
!!$
!!$ end do
!!$
!!$ return
!!$ end subroutine DensityProfile
!!$
!!$!
!!$
!!$ subroutine PrintDensity(dens)
!!$
!!$ implicit none
!!$
!!$ real (kind=8) :: x,y
!!$ integer (kind=4) :: i,j
!!$
!!$ real (kind=8), dimension (Nbin,Nbin) :: dens
!!$
!!$ do j=1,Nbin
!!$ y = -0.5*rcut+j*rbin
!!$ do i=1,Nbin
!!$ x = -0.5*rcut+i*rbin
!!$ write (101,*) x,y,dens(i,j)/rbin**2
!!$ end do
!!$ write (101,*)
!!$ end do
!!$
!!$ return
!!$ end subroutine PrintDensity
!-----------------------------------------------------------------------
subroutine NormalizeGr(density,ngr,gr)
implicit none
real (kind=8) :: r8_gamma
real (kind=8) :: density
real (kind=8) :: nid,r,norm
real (kind=8) :: k_n
integer (kind=4) :: ibin
integer (kind=4) :: ngr
real (kind=8),dimension (Nbin) :: gr
k_n = pi**(0.5d0*dim)/r8_gamma(0.5d0*dim+1.d0)
norm = real(Np)*real(ngr)
do ibin=1,Nbin
r = (real(ibin)-0.5d0)*rbin
nid = density*k_n*((r+0.5d0*rbin)**dim-(r-0.5d0*rbin)**dim)
gr(ibin) = gr(ibin)/(nid*norm)
end do
return
end subroutine NormalizeGr
!-----------------------------------------------------------------------
subroutine NormalizeSk(Nk,ngr,Sk)
implicit none
real (kind=8) :: norm
integer (kind=4) :: k,ik,Nk
integer (kind=4) :: ngr
real (kind=8), dimension (dim,Nk) :: Sk
norm = real(Np)*real(ngr)
do ik=1,Nk
do k=1,dim
Sk(k,ik) = Sk(k,ik)/norm
end do
end do
return
end subroutine NormalizeSk
!-----------------------------------------------------------------------
subroutine NormalizeNr(density,zconf,Nobdm,nrho)
implicit none
real (kind=8) :: r8_gamma
real (kind=8) :: density
real (kind=8) :: k_n,nid
real (kind=8) :: zconf
real (kind=8) :: r
integer (kind=4) :: ibin,m
integer (kind=4) :: Nobdm
real (kind=8), dimension (0:Npw,Nbin) :: nrho
k_n = pi**(0.5d0*dim)/r8_gamma(0.5d0*dim+1.d0)
do ibin=1,Nbin
r = (real(ibin)-0.5d0)*rbin
nid = density*k_n*((r+0.5d0*rbin)**dim-(r-0.5d0*rbin)**dim)
do m=0,Npw
nrho(m,ibin) = nrho(m,ibin)/(CWorm*nid*zconf*real(Nobdm))
end do
end do
return
end subroutine NormalizeNr
!-----------------------------------------------------------------------
subroutine AccumGr(gr,AvGr,AvGr2)
implicit none
integer (kind=4) :: j
real (kind=8),dimension(Nbin) :: gr,AvGr,AvGr2
do j=1,Nbin
AvGr(j) = AvGr(j)+gr(j)
AvGr2(j) = AvGr2(j)+gr(j)*gr(j)
end do
return
end subroutine AccumGr
!-----------------------------------------------------------------------
subroutine AccumSk(Nk,Sk,AvSk,AvSk2)
implicit none
integer (kind=4) :: j,k,Nk
real (kind=8),dimension(dim,Nk) :: Sk,AvSk,AvSk2
do j=1,Nk
do k=1,dim
AvSk(k,j) = AvSk(k,j)+Sk(k,j)
AvSk2(k,j) = AvSk2(k,j)+Sk(k,j)*Sk(k,j)
end do
end do
return
end subroutine AccumSk
!-----------------------------------------------------------------------
subroutine AccumNr(nrho,AvNr,AvNr2)
implicit none
integer (kind=4) :: j,m
real (kind=8),dimension(0:Npw,Nbin) :: nrho,AvNr,AvNr2
do j=1,Nbin
do m=0,Npw
AvNr(m,j) = AvNr(m,j)+nrho(m,j)
AvNr2(m,j) = AvNr2(m,j)+nrho(m,j)*nrho(m,j)
end do
end do
return
end subroutine AccumNr
!-----------------------------------------------------------------------
subroutine NormAvGr(Nitem,AvGr,AvGr2,VarGr)
implicit none
real (kind=8) :: r
integer (kind=4) :: Nitem,j
real (kind=8),dimension(Nbin) :: AvGr,AvGr2,VarGr
open (unit=11,file='gr_vpi.out')
do j=1,Nbin
r = (real(j)-0.5d0)*rbin
AvGr(j) = AvGr(j)/real(Nitem)
AvGr2(j) = AvGr2(j)/real(Nitem)
VarGr(j) = Var(Nitem,AvGr(j),AvGr2(j))
write (11,'(20g20.10e3)') r,AvGr(j),VarGr(j)
end do
close (unit=11)
return
end subroutine NormAvGr
!-----------------------------------------------------------------------
subroutine NormAvSk(Nitem,Nk,AvSk,AvSk2,VarSk)
implicit none
integer (kind=4) :: Nitem,Nk,j,k
real (kind=8),dimension(dim,Nk) :: AvSk,AvSk2,VarSk
open (unit=11,file='sk_vpi.out')
do j=1,Nk
do k=1,dim
AvSk(k,j) = AvSk(k,j)/real(Nitem)
AvSk2(k,j) = AvSk2(k,j)/real(Nitem)
VarSk(k,j) = Var(Nitem,AvSk(k,j),AvSk2(k,j))
end do
write (11,'(20g20.10e3)') (j*qbin(k),AvSk(k,j),VarSk(k,j),k=1,dim)
end do
close (unit=11)
return
end subroutine NormAvSk
!-----------------------------------------------------------------------
subroutine NormAvNr(Nitem,AvNr,AvNr2,VarNr)
implicit none
real (kind=8) :: r
integer (kind=4) :: Nitem,j,m
real (kind=8),dimension(0:Npw,Nbin) :: AvNr,AvNr2,VarNr
open (unit=11,file='nr_vpi.out')
do j=1,Nbin
r = (real(j)-0.5d0)*rbin
do m=0,Npw
AvNr(m,j) = AvNr(m,j)/real(Nitem)
AvNr2(m,j) = AvNr2(m,j)/real(Nitem)
VarNr(m,j) = Var(Nitem,AvNr(m,j),AvNr2(m,j))
end do
write (11,'(20g20.10e3)') r,(AvNr(m,j),VarNr(m,j),m=0,Npw)
end do
close (unit=11)
return
end subroutine NormAvNr
!-----------------------------------------------------------------------
!
! This subroutine simply accumulates the value of the energies at each
! step/block in order to obtain average values.
!
!-----------------------------------------------------------------------
subroutine Accumulate(E,Ec,Ep,SumE,SumEc,SumEp)
implicit none
real (kind=8) :: E,SumE
real (kind=8) :: Ec,SumEc
real (kind=8) :: Ep,SumEp
SumE = SumE+E
SumEc = SumEc+Ec
SumEp = SumEp+Ep
return
end subroutine Accumulate
!-----------------------------------------------------------------------
!
! This subroutine evaluates the average of a magnitude along a step,
! block or the whole simulation.
!
!-----------------------------------------------------------------------
subroutine NormalizeAv(Nitem,SumE,SumEc,SumEp)
implicit none
real (kind=8) :: SumE,SumEc,SumEp
integer (kind=4) :: Nitem
SumE = SumE/real(Nitem)
SumEc = SumEc/real(Nitem)
SumEp = SumEp/real(Nitem)
return
end subroutine NormalizeAv
!-----------------------------------------------------------------------
!
! This function simply evaluates the variance of a variable.
!
!-----------------------------------------------------------------------
function Var(Nitem,Sum,Sum2)
implicit none
real (kind=8) :: Var
real (kind=8) :: Sum,Sum2
integer (kind=4) :: Nitem
Var = sqrt((Sum2-Sum*Sum)/real(Nitem))
return
end function Var
!-----------------------------------------------------------------------
end module sample_mod