-
Notifications
You must be signed in to change notification settings - Fork 2
/
oz_mod.f90
1892 lines (1547 loc) · 61.8 KB
/
oz_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
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
! This file is part of SunlightDPD - a home for open source software
! related to the dissipative particle dynamics (DPD) simulation
! method.
! Based on an original code copyright (c) 2007 Lucian Anton.
! Modifications copyright (c) 2008, 2009 Andrey Vlasov. Additional
! modifications copyright (c) 2009-2018 Unilever UK Central Resources
! Ltd (Registered in England & Wales, Company No 29140; Registered
! Office: Unilever House, Blackfriars, London, EC4P 4BQ, UK). Later
! modifications copyright (c) 2020-2024 Patrick B Warren
! <patrick.warren@stfc.ac.uk> and STFC.
! SunlightDPD is free software: you can redistribute it and/or
! modify it under the terms of the GNU General Public License as
! published by the Free Software Foundation, either version 3 of the
! License, or (at your option) any later version.
! SunlightDPD is distributed in the hope that it will be useful, but
! WITHOUT ANY WARRANTY; without even the implied warranty of
! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
! General Public License for more details.
! You should have received a copy of the GNU General Public License
! along with SunlightDPD. If not, see <http://www.gnu.org/licenses/>.
! Full documentation is found in the accompanying LaTeX document.
module wizard
implicit none
include "fftw3.f"
integer, parameter :: dp = kind(1.0d0)
integer, parameter :: & ! Enumeration of current closures
& NO_CLOSURE = 0, &
& HNC_CLOSURE = 1, &
& RPA_CLOSURE = 2, &
& MSA_CLOSURE = 3, &
& EXP_CLOSURE = 4
integer, parameter :: & ! Enumeration of current potential models
& NO_MODEL_TYPE = 0, &
& DPD_GAUSSIAN_CHARGES = 1, &
& DPD_BESSEL_CHARGES = 2, &
& DPD_LINEAR_CHARGES = 3, &
& DPD_SLATER_APPROX_CHARGES = 4, &
& DPD_SLATER_EXACT_CHARGES = 5, &
& URPM_WITHOUT_USHORT = 6, &
& URPM_WITH_USHORT = 7, &
& RPM_WITHOUT_USHORT = 8, &
& RPM_WITH_USHORT = 9, &
& HARD_SPHERES = 10
integer, parameter :: & ! Enumeration of error codes
& NO_ERROR = 0, &
& CONVERGENCE_ERROR = 1, &
& AXEQB_ERROR = 2, &
& DSYSV_ERROR = 3, &
& DGEEV_ERROR = 4, &
& MISMATCH_ERROR = 5
character (len=5) :: version = 'v1.13' ! The current version
character (len=3) :: closure_name = '' ! TLA for the last-used closure
character (len=32) :: model_name = '' ! Model name
character (len=48) :: error_msg = '' ! Error message
logical :: suppress_msgs = .false. ! In case of severe numerical instability
logical :: silent = .false. ! Prevent printing warning/error messages
logical :: verbose = .false. ! Print solver diagnostics
logical :: cold_start = .true. ! Force a cold start
logical :: auto_fns = .true. ! Calculate stuff at end
real(kind=dp), parameter :: &
& pi = 4.0_dp * atan(1.0_dp), &
& twopi = 2.0_dp * pi, &
& fourpi = 4.0_dp * pi, &
& rootpi = sqrt(pi)
integer :: &
& start_type = 3, & ! how to initialise in a cold start
& model_type = 0, & ! which potential was (last) chosen
& closure_type = 0, & ! last-used closure
& istep, & ! current position in iterative solver
& ng = 4096, & ! grid size
& ncomp = 1, & ! number of chemical components
& nfnc = 0, & ! number of functions, = ncomp (ncomp + 1) / 2
& nps = 6, & ! number of previous states used in Ng method
& npic = 6, & ! number of Picard steps
& maxsteps = 100, & ! max number of steps to take for convergence
& nsteps, & ! actual number of steps taken
& return_code = 0 ! error code (see above)
integer*8 :: plan ! FFTW plan for fast discrete sine transforms
real(kind=dp) :: &
& deltar = 0.01_dp, & ! real space grid spacing
& deltak, & ! reciprocal space grid spacing (computed)
& error, & ! difference between current and previous solns
& alpha = 0.2_dp, & ! Picard method, fraction of new solution
& tol = 1.0E-12_dp, & ! Error tolerance for claiming convergence
& rc = 1.0_dp, & ! short-range DPD repulsion range
& lb = 0.0_dp, & ! long-range Coulomb coupling length
& sigma = 1.0_dp, & ! see below (*)
& sigmap = 1.0_dp, & ! +- long-range Coulomb smearing length (URPM)
& kappa = -1.0_dp, & ! +- long-range Coulomb smoothing parameter (RPM)
& rgroot = 1.0_dp, & ! linear charge smearing range (Groot)
& lbda = 1.0_dp, & ! Slater charge smearing range (exact)
& beta = 1.0_dp, & ! Slater charge smearing range (approx)
& cf_mf, cf_xc, & ! the virial route pressure contributions ..
& cf_gc, pex, press, & ! .. and the virial route pressure
& comp, comp_xc, & ! compressibility, and excess
& aex_rl, aex_rs, & ! contributions to HNC excess free energy density
& aex_kl, aex_kd, & ! --- '' ---
& aex_ks, aex, & ! -- '' --, and final excess free energy density
& deficit, & ! excess free energy deficit (see docs)
& uex, uex_mf, uex_xc ! energy density total and contributions
! (*) sigma is used both for the long-range Coulomb smearing length
! for the soft potentials, and for hard core diameter for RPM models
real(kind=dp), allocatable :: &
& rho(:), & ! density array
& z(:), & ! valence array
& arep(:, :), & ! repulsion amplitude array
& dd(:), & ! hard core diameter array
& diam(:, :), & ! ditto, pairwise
& tp(:), & ! mean field pressure contribution
& tu(:), & ! mean field energy contribution
& tl(:), & ! mean field long range potential
& muex(:), & ! chemical potential array
& c(:, :, :), & ! direct correlation functions (dcfs)
& e(:, :, :), & ! indirect correlation functions (icfs)
& h0(:, :), & ! reference total correlation functions
& hr(:, :, :), & ! current total correlation functions
& hc(:, :), & ! current contact values
& ck(:, :), & ! transform of dcfs
& ek(:, :), & ! transform of icfs
& sk(:, :, :), & ! partial structure factors
& ushort(:, :), & ! short range potential in real space
& expnegus(:, :), & ! exp(-ushort) (includes hard cores)
& dushort(:, :), & ! derivative of the same
& ulong(:, :), & ! long range potential in real space
& dulong(:, :), & ! derivative of the same
& ulongk(:, :), & ! long range potential in reciprocal space
& u0(:), & ! r = 0 of long range potential (diagonal)
& r(:), k(:), & ! r and k grids
& fftwx(:), fftwy(:) ! arrays for fast discrete sine transform
contains
subroutine initialise
implicit none
integer i
nfnc = ncomp * (ncomp + 1) / 2
allocate(rho(ncomp))
allocate(z(ncomp))
allocate(muex(ncomp))
allocate(arep(ncomp, ncomp))
allocate(diam(ncomp, ncomp))
allocate(dd(nfnc))
allocate(tp(nfnc))
allocate(tu(nfnc))
allocate(tl(nfnc))
allocate(c(ng-1, nfnc, nps))
allocate(e(ng-1, nfnc, nps))
allocate(h0(ng-1, nfnc))
allocate(hr(ng-1, ncomp, ncomp))
allocate(hc(ncomp, ncomp))
allocate(ck(ng-1, nfnc))
allocate(ek(ng-1, nfnc))
allocate(sk(ng-1, ncomp, ncomp))
allocate(ushort(ng-1, nfnc))
allocate(expnegus(ng-1, nfnc))
allocate(dushort(ng-1, nfnc))
allocate(ulong(ng-1, nfnc))
allocate(dulong(ng-1, nfnc))
allocate(ulongk(ng-1, nfnc))
allocate(u0(ncomp))
allocate(r(ng-1))
allocate(k(ng-1))
allocate(fftwx(ng-1))
allocate(fftwy(ng-1))
! Default values
rho = 0.0_dp
arep = 0.0_dp
diam = 0.0_dp
z = 0.0_dp
dd = 0.0_dp
h0 = 0.0_dp
! Make grids
deltak = pi / (real(ng, kind=dp) * deltar)
r = (/ (i*deltar, i=1, ng-1) /)
k = (/ (i*deltak, i=1, ng-1) /)
! Make the FFTW plan
call dfftw_plan_r2r_1d(plan, ng-1, fftwx, fftwy, &
& FFTW_RODFT00, FFTW_ESTIMATE)
end subroutine initialise
! Write out parameters for system and potentials.
subroutine write_params
implicit none
integer :: i
print *, '====================================================='
print *, 'GRID DETAILS'
print *, ' ng = ', ng, ' ncomp = ', ncomp, ' nfnc = ', nfnc, ' nps = ', nps
print *, ' deltar = ', deltar, ' deltak = ', deltak
print *, ' deltar*deltak*ng/pi = ', deltar*deltak/pi*dble(ng)
print *, ' r(ng-1) = ', r(ng-1), ' k(ng-1) = ', k(ng-1)
print *, 'POTENTIAL DETAILS :: model type =', model_type, model_name
if (model_type.eq.NO_MODEL_TYPE) then
print *, 'No potential has been selected / potential set externally'
else if (model_type.le.DPD_SLATER_EXACT_CHARGES) then
print *, 'DPD potential was selected, matrix A = '
do i = 1, ncomp
print *, ' ', arep(i, :)
end do
print *, ' valencies, z = ', z
print *, ' rc = ', rc, ' lb = ', lb
if (model_type.eq.DPD_GAUSSIAN_CHARGES) then
print *, ' Gaussian smearing, sigma = ', sigma
end if
if (model_type.eq.DPD_BESSEL_CHARGES) then
print *, ' Bessel smearing, sigma = ', sigma
end if
if (model_type.eq.DPD_LINEAR_CHARGES) then
print *, ' linear (Groot) smearing, R = ', rgroot
print *, ' equivalent Gaussian sigma = ', &
& sqrt(2.0_dp/15.0_dp) * rgroot
end if
if (model_type.eq.DPD_SLATER_APPROX_CHARGES) then
print *, ' Slater smearing (approx), beta = ', beta
print *, ' 1 / beta = ', 1.0_dp / beta, ', &
& 5 / (8 beta) = ', 0.625_dp / beta
print *, ' equivalent Gaussian sigma = ', sqrt(0.5_dp) / beta
end if
if (model_type.eq.DPD_SLATER_EXACT_CHARGES) then
print *, ' Slater smearing (exact), lambda = ', lbda
print *, ' equivalent Gaussian sigma = ', lbda
end if
else if (model_type.le.URPM_WITH_USHORT) then
if (model_type.eq.URPM_WITHOUT_USHORT) then
print *, 'URPM potential was selected with ushort unused'
else
print *, 'URPM potential was selected with ushort used'
end if
print *, ' lb = ', lb, ' sigma = ', sigma, ' sigmap = ', sigmap
else if (model_type.le.RPM_WITH_USHORT) then
if (model_type.eq.RPM_WITHOUT_USHORT) then
print *, 'RPM potential was selected with ushort unused'
else
print *, 'RPM potential was selected with ushort used'
end if
if (kappa.lt.0.0_dp) then
print *, ' lb = ', lb, ' sigma = ', sigma, ' kappa -> infinity'
else
print *, ' lb = ', lb, ' sigma = ', sigma, ' kappa = ', kappa
end if
print *, ' valencies, z = ', z
print *, 'hard core diameter array'
do i = 1, ncomp
print *, ' ', diam(i, :)
end do
print *, ' hard sphere diameters, dd = ', dd
else if (model_type.eq.HARD_SPHERES) then
print *, 'Hard sphere potential with sigma = ', sigma
print *, 'hard core diameter array'
do i = 1, ncomp
print *, ' ', diam(i, :)
end do
print *, ' hard sphere diameters, dd = ', dd
else
print *, 'Undefined potential'
end if
print *, 'SYSTEM DETAILS'
print *, ' rho = ', rho
if (sum(rho).gt.0.0_dp) print *, ' x = ', rho(:) / sum(rho)
print *, ' sum(rho) = ', sum(rho)
print *, ' sum(rho*z) = ', sum(rho(:)*z(:))
print *, '====================================================='
end subroutine write_params
! Build the potential arrays, with parameters rc and arep(:,:) for the
! short-range DPD repulsion, and lb, sigma and z(:) for the long-range
! Gaussian-smeared Coulomb part. A factor beta = 1/kT is implicit in
! these definitions. The parameter charge_type is (default) 1 =
! Gaussian, 2 = Bessel, 3 = Groot (linear), 4 = Slater (approx), 5 =
! Slater (exact). But use the defined integer constants for these!
subroutine dpd_potential(charge_type)
implicit none
integer, intent(in), optional :: charge_type
integer :: ctype = DPD_GAUSSIAN_CHARGES
integer :: i, j, ij, irc
real(kind=dp) :: ulong0, aa(nfnc), zz(nfnc)
if (present(charge_type) .and. charge_type.gt.0) ctype = charge_type
model_type = ctype
! Sort out some recoded potential parameters. Also force the
! amplitude matrix to be symmetric, set by upper triangular
! entries. This enforces the rule in the documentation and also
! simplifies the printing.
do j = 1, ncomp
do i = 1, j
ij = i + j*(j-1)/2
aa(ij) = arep(i, j)
zz(ij) = z(i) * z(j)
if (i.lt.j) arep(j, i) = arep(i, j)
end do
end do
irc = nint(rc/deltar)
! Leave out the amplitude, then the function can be re-used
! (see below)
ushort(:,nfnc) = 0.0_dp
ushort(1:irc,nfnc) = 0.5_dp * (1.0_dp - r(1:irc)/rc)**2
dushort(:,nfnc) = 0.0_dp
dushort(1:irc,nfnc) = - (1.0_dp - r(1:irc)/rc) / rc
! Gaussian charges
if (ctype .eq. DPD_GAUSSIAN_CHARGES) then
ulong(:,nfnc) = (lb / r) * erf(0.5_dp*r/sigma)
ulong0 = lb / (sigma * rootpi)
ulongk(:,nfnc) = (fourpi * lb / k**2) * exp(-k**2*sigma**2)
dulong(:,nfnc) = lb * exp(-0.25_dp*r**2/sigma**2) &
& / (rootpi * r * sigma) &
& - lb * erf(0.5_dp*r/sigma) / r**2
model_name = "DPD_with_Gaussian_charges"
end if
! Bessel charges
if (ctype .eq. DPD_BESSEL_CHARGES) then
ulong(:,nfnc) = (lb / r) * (1.0_dp - exp(-r/sigma))
ulong0 = lb / sigma
ulongk(:,nfnc) = (fourpi * lb / k**2) &
& * 1.0_dp / (1.0_dp + k**2*sigma**2)
dulong(:,nfnc) = - (lb / r**2) * (1.0_dp - exp(-r/sigma) &
& * (1.0_dp + r / sigma))
model_name = "DPD_with_Bessel_charges"
end if
! Linear charge smearing as in Groot [JCP v118, 11265 (2003)].
! Note we do not give the real space part here hence the
! thermodynamic calculations will be wrong. This could be fixed
! up by doing a numerical FFT of the potential. Best would be to
! separate off the long range 4 pi lb / k^2 first. TO BE DONE!
if (ctype .eq. DPD_LINEAR_CHARGES) then
ulong(:,nfnc) = 0.0_dp
ulong0 = 0.0_dp
ulongk(:,nfnc) = (fourpi * lb / k**2) * 144.0_dp &
& * (2.0_dp - 2.0_dp*cos(k*rgroot) &
& - k*rgroot*sin(k*rgroot))**2 &
& / (k**8 * rgroot**8)
dulong(:,nfnc) = 0.0_dp
model_name = "DPD_with_linear_charges"
end if
! Slater charge smearing as in Gonzales-Melchor et al, [JCP v125,
! 224107 (2006)] but here with exact expression for interaction.
if (ctype .eq. DPD_SLATER_EXACT_CHARGES) then
ulong(:,nfnc) = (lb / r) * (1.0_dp - exp(-2.0_dp*r/lbda) &
& * (1.0_dp + 1.375_dp*r/lbda + 0.75_dp*r**2/lbda**2 &
& + r**3/(6.0_dp*lbda**3)) )
ulong0 = 0.625_dp * lb / lbda
ulongk(:,nfnc) = (fourpi * lb / k**2) * &
& 1.0_dp / (1.0_dp + k**2*lbda**2/4.0_dp)**4
dulong(:,nfnc) = - (lb / r**2) * (1.0_dp - exp(-2.0_dp*r/sigma) &
& * (1.0_dp + 2.0_dp*r/lbda + 2.0_dp*r**2/lbda**2 &
& + 7.0_dp*r**3/(6.0_dp*lbda**3) + r**4/(3.0_dp*lbda**4)) )
model_name = "DPD_with_Slater_(exact)_charges"
end if
! Slater charge smearing as in Gonzales-Melchor et al, [JCP v125,
! 224107 (2006)] with original approximate expression.
if (ctype .eq. DPD_SLATER_APPROX_CHARGES) then
ulong(:,nfnc) = (lb / r) * (1.0_dp - exp(-2*beta*r) * &
& (1.0_dp + beta*r) )
ulong0 = beta * lb
ulongk(:,nfnc) = (fourpi * lb / k**2) * &
& 1.0_dp / (1.0_dp + k**2/(4.0_dp*beta**2))**2
dulong(:,nfnc) = - (lb / r**2) * (1.0_dp - exp(-2.0_dp*beta*r) &
& * (1.0_dp + 2.0_dp*beta*r*(1.0_dp + beta*r)) )
model_name = "DPD_with_Slater_(approx)_charges"
end if
! Generate the pair potentials by walking through the functions.
! In the final step we correctly normalise the final function.
do i = 1, nfnc
ushort(:,i) = aa(i) * ushort(:,nfnc)
dushort(:,i) = aa(i) * dushort(:,nfnc)
ulong(:,i) = zz(i) * ulong(:,nfnc)
ulongk(:,i) = zz(i) * ulongk(:,nfnc)
dulong(:,i) = zz(i) * dulong(:,nfnc)
end do
! These individual species-pair contributions to the mean field
! compressibility factor and the mean-field internal energy per
! particle can be calculated analytically for the DPD potential.
tp = pi * rc**3 * aa / 30.0_dp
tu = tp
tl = 0.0_dp
! The r = 0 diagonal parts of the long range potential
u0(:) = z(:)**2 * ulong0
! Generate auxiliary function
expnegus = exp(-ushort)
end subroutine dpd_potential
! Build the potential arrays for the softened URPM (Gaussian charges),
! with parameters lb, sigma and sigmap. This expects ncomp = 2, and
! will set z(1) = 1, z(2) = -1. The logical parameter controls
! whether ushort is used or not.
subroutine urpm_potential(use_ushort)
implicit none
logical, intent(in), optional :: use_ushort
logical :: uuflag = .false.
if (present(use_ushort)) uuflag = use_ushort
if (ncomp.ne.2) then
return_code = MISMATCH_ERROR
error_msg = 'mismatch ncomp <> 2 in urpm_potential'
if (.not.silent) print *, '** error: ', error_msg
return
end if
z(1) = 1.0_dp
z(2) = -1.0_dp
ulong(:,1) = lb * erf(0.5_dp*r/sigma) / r
ulongk(:,1) = fourpi * lb * exp(-k**2*sigma**2) / k**2
dulong(:,1) = lb * exp(-0.25_dp*r**2/sigma**2) / (rootpi * r * sigma) &
& - lb * erf(0.5_dp*r/sigma) / r**2
ulong(:,2) = - lb * erf(0.5_dp*r/sigmap) / r
ulongk(:,2) = - fourpi * lb * exp(-k**2*sigmap**2) / k**2
dulong(:,2) = - lb * exp(-0.25_dp*r**2/sigmap**2) / (rootpi * r * sigmap) &
& + lb * erf(0.5_dp*r/sigmap) / r**2
ulong(:,3) = ulong(:,1)
ulongk(:,3) = ulongk(:,1)
dulong(:,3) = dulong(:,1)
ushort(:,:) = 0.0_dp
dushort(:,:) = 0.0_dp
if (uuflag) then
ushort(:,2) = ulong(:,2) + ulong(:,1)
dushort(:,2) = dulong(:,2) + dulong(:,1)
ulong(:,2) = - ulong(:,1)
ulongk(:,2) = - ulongk(:,1)
dulong(:,2) = - dulong(:,1)
end if
! These individual species-pair contributions to the mean field
! compressibility factor and the mean-field internal energy per
! particle can be calculated analytically for the URPM potential.
! These are the same whether using ushort or not, as they are
! defined in terms of the total potential.
tp(1) = 0.0_dp
tp(2) = twopi*lb*(sigmap**2 - sigma**2)
tp(3) = 0.0_dp
tu = tp
! If not using ushort, we are off the symmetry point condition and
! the contribution of the long range part should be incorporated
! into the compressibility and chemical potential expressions.
if (.not.uuflag) then
tl = 2.0_dp*tp
else
tl = 0.0_dp
end if
! The r = 0 diagonal parts of the long range potential
u0(:) = z(:)**2 * lb / (sigma * rootpi)
! Generate auxiliary function
expnegus = exp(-ushort)
if (.not.uuflag) then
model_type = URPM_WITHOUT_USHORT
model_name = 'URPM_without_U_short'
else
model_type = URPM_WITH_USHORT
model_name = 'URPM_with_U_short'
end if
end subroutine urpm_potential
! Build the potential arrays for the softened RPM (charged hard
! spheres) with parameters lb, sigma and kappa. This expects ncomp =
! 2 or 3, and will set z(1) = 1, z(2) = -1, and z(3) = 0 if required,
! and hard core diameters to sigma. The parameter (0 or 1) controls
! whether ushort is used or not. Using kappa < 0 implies the pure RPM
! case (kappa -> infinity). The case ncomp = 3 corresponds to the RPM
! in the presence of a neutral hard sphere solvent.
subroutine rpm_potential(use_ushort)
implicit none
integer :: i, j, ij, irc
logical, intent(in), optional :: use_ushort
logical :: uuflag = .false.
if (present(use_ushort)) uuflag = use_ushort
if (ncomp.eq.1 .or. ncomp.gt.3) then
return_code = MISMATCH_ERROR
error_msg = 'mismatch ncomp <> 2 or 3 in rpm_potential'
if (.not.silent) print *, '** error: ', error_msg
return
end if
z(1) = 1.0_dp
z(2) = -1.0_dp
if (ncomp.eq.3) z(3) = 0.0_dp
! Set the diameters from an array, or sigma if none are preset.
! The remaining logic means sigma, if not set, is set to the
! smallest diameter
if (any(diam.gt.0.0_dp)) then
do j = 1, ncomp
do i = 1, j
ij = i + j*(j-1)/2
dd(ij) = diam(i, j)
end do
end do
if (sigma.eq.0.0_dp) sigma = minval(dd)
else
dd = sigma
end if
! Zero everything to start with, as for hard spheres
ulong = 0.0_dp
ulongk = 0.0_dp
dulong = 0.0_dp
ushort = 0.0_dp
dushort = 0.0_dp
irc = nint(sigma / deltar)
ulong(:,1) = lb / r
ulong(1:irc,1) = lb / sigma
ulongk(:,1) = fourpi * lb * sin(k*sigma) / (sigma * k**3)
dulong(:,1) = - lb / r**2
dulong(1:irc,1) = 0.0_dp
if (kappa.gt.0.0_dp) then
ulong(:,2) = - lb * erf(kappa*r) / r
ulongk(:,2) = - fourpi * lb * exp(-k**2/(4.0_dp*kappa**2)) / k**2
dulong(:,2) = - 2.0_dp*kappa*lb * exp(-kappa**2*r**2) / (rootpi * r) &
& + lb * erf(kappa*r) / r**2
else
ulong(:,2) = - ulong(:,1)
ulongk(:,2) = - ulongk(:,1)
dulong(:,2) = - dulong(:,1)
end if
ulong(:,3) = ulong(:,1)
ulongk(:,3) = ulongk(:,1)
dulong(:,3) = dulong(:,1)
ushort(:,:) = 0.0_dp
dushort(:,:) = 0.0_dp
if (uuflag) then
ushort(:,2) = ulong(:,2) + ulong(:,1)
dushort(:,2) = dulong(:,2) + dulong(:,1)
ulong(:,2) = - ulong(:,1)
ulongk(:,2) = - ulongk(:,1)
dulong(:,2) = - dulong(:,1)
end if
! Generate auxiliary function - this is where HS diam implemented
expnegus = exp(-ushort)
do ij = 1, nfnc
irc = nint(dd(ij) / deltar)
expnegus(1:irc, ij) = 0.0_dp
end do
! These are the analytic contributions to the thermodynamics.
tp = 0.0_dp
tu = 0.0_dp
tl = 0.0_dp
if (kappa.gt.0) then
tp(2) = pi*lb * ( sigma * exp(-kappa**2*sigma**2) / (kappa*rootpi) &
& + (1/(2.0_dp*kappa**2) - sigma**2/3.0_dp) * erfc(kappa*sigma) )
tu(2) = pi*lb * ( sigma * exp(-kappa**2*sigma**2) / (kappa*rootpi) &
& + (1/(2.0_dp*kappa**2) - sigma**2) * erfc(kappa*sigma) )
if (.not.uuflag) then ! off SYM condition;
! missing last term added in v1.11
tl(2) = pi*lb / kappa**2 - 2.0_dp*pi*lb*sigma**2 / 3.0_dp
end if
end if
! The r = 0 diagonal parts of the long range potential
u0(:) = z(:)**2 * lb / sigma
if (.not.uuflag) then
model_type = RPM_WITHOUT_USHORT
model_name = 'RPM_without_U_short'
else
model_type = RPM_WITH_USHORT
model_name = 'RPM_with_U_short'
end if
end subroutine rpm_potential
! Build the potential arrays for hard spheres with diameter sigma.
! This works for arbitrary numbers of components.
subroutine hs_potential
implicit none
integer :: i, j, ij, irc
! Set the hard core diameters from an array, or sigma if none are pre-set
if (any(diam.gt.0.0_dp)) then
do j = 1, ncomp
do i = 1, j
ij = i + j*(j-1)/2
dd(ij) = diam(i, j)
end do
end do
else
dd = sigma
end if
ulong = 0.0_dp
ulongk = 0.0_dp
dulong = 0.0_dp
ushort = 0.0_dp
dushort = 0.0_dp
! This is where the hard sphere diameters enter
expnegus = 1.0_dp
do ij = 1, nfnc
irc = nint(dd(ij) / deltar)
expnegus(1:irc, ij) = 0.0_dp
end do
! These are the analytic contributions to the thermodynamics.
tp = 0.0_dp
tu = 0.0_dp
tl = 0.0_dp
u0 = 0.0_dp
model_type = HARD_SPHERES
model_name = 'hard_spheres'
end subroutine hs_potential
! The next routine solves the Ornstein-Zernicke equation to determine
! e = h - c, given c. We re-partition the long range part of the
! potential so that the routine actually calculates c' = c + Ulong
! and e' = e - Ulong. This is because the Fourier transform of
! Ulong can be computed in closed form. Note h = e + c = e' + c'.
subroutine oz_solve
implicit none
integer :: i1, i, j, ij, ik, irc
integer :: perm(ncomp)
real(kind=dp) :: &
& a(ncomp, ncomp), b(ncomp, ncomp), &
& cmat(ncomp, ncomp), umat(ncomp, ncomp), rhomat(ncomp, ncomp), &
& m0(ncomp, ncomp), unita(ncomp, ncomp)
i1 = mod(istep-1, nps) + 1
! Forward transform the real space functions c, to the reciprocal
! space functions ck.
do i=1, nfnc
fftwx(1:ng-1) = r(1:ng-1) * c(1:ng-1, i, i1)
call dfftw_execute(plan)
ck(1:ng-1, i) = (twopi * deltar) * fftwy(1:ng-1) / k(1:ng-1)
end do
if (ncomp .eq. 1) then
! In the one component case OZ inversion is straightforward.
! Note the implicit indexing on wavevector k.
ek(:, 1) = ( ck(:, 1) - ulongk(:, 1) ) &
& / ( 1.0_dp - rho(1) * (ck(:, 1) - ulongk(:, 1)) ) &
& - ck(:, 1)
else ! Multicomponent OZ inversion
! First set up a unit matrix, and the diagonal R matrix -- see
! the documentation for the math here.
rhomat = 0.0_dp
unita = 0.0_dp
do i = 1, ncomp
rhomat(i,i) = rho(i)
unita(i,i) = 1.0_dp
end do
! Do the matrix calculations for each wavevector k.
do ik = 1, ng-1
! Unpack the reciprocal space functions into matrices.
do j = 1, ncomp
do i = 1, j
ij = i + j*(j-1)/2
cmat(i, j) = ck(ik, ij)
umat(i, j) = ulongk(ik, ij)
if (i.lt.j) then
cmat(j, i) = cmat(i, j)
umat(j, i) = umat(i, j)
end if
end do
end do
! The following matrices are constructed:
! M0 = (C - beta UL) . R
! A = I - (C - beta UL) . R
! B = (C - beta UL) . R . C - beta UL
! (note that beta = 1/kT = 1 is not written explicitly)
m0 = matmul(cmat - umat, rhomat)
a = unita - m0
b = matmul(m0, cmat) - umat
! Solve the equation A.X = B so that
! X = [I - (C - beta U) . R]^(-1) . [(C - beta U) . R . C - beta U]
call axeqb_reduce(a, ncomp, b, ncomp, perm, irc)
if (irc.gt.0) then
return_code = AXEQB_ERROR
error_msg = 'axeqb encountered singular problem in oz_solve'
if (.not.suppress_msgs .and. .not.silent) then
print *, '** error: ', error_msg
print *, '** further messages of this kind will be suppressed'
suppress_msgs = .true.
end if
return
end if
! Now X(I, :) = B(PERM(I), :) is the new estimate for the
! reciprocal space functions ek. They are built from the
! upper triangle of the matrix.
do j = 1, ncomp
do i = 1, j
ij = i + j*(j-1)/2
ek(ik, ij) = b(perm(i), j)
end do
end do
end do ! loop over k vectors
end if ! select single component or multicomponent case
! Do the Fourier back transforms
do i = 1, nfnc
fftwx(1:ng-1) = k(1:ng-1) * ek(1:ng-1, i)
call dfftw_execute(plan)
e(1:ng-1, i, i1) = (deltak / twopi**2) * fftwy(1:ng-1) / r(1:ng-1)
end do
end subroutine oz_solve
! This routine solves an alternate version of the Ornstein-Zernicke
! equation to determine c and e from the reference h0. In practice as
! always we actually calculate c' = c + Ulong and e' = e - Ulong.
! Note h = e + c = e' + c'. The result is saved to position 1 in the
! history trajectory. Note that the offset ulongk in
subroutine oz_solve2
implicit none
integer :: i, j, ij, ik, irc
integer :: perm(ncomp)
real(kind=dp) :: &
& a(ncomp, ncomp), b(ncomp, ncomp), &
& hmat(ncomp, ncomp), rhomat(ncomp, ncomp), &
& unita(ncomp, ncomp), hk(ng-1, nfnc)
do i=1, nfnc
fftwx(1:ng-1) = r(1:ng-1) * h0(1:ng-1, i)
call dfftw_execute(plan)
hk(1:ng-1, i) = (twopi * deltar) * fftwy(1:ng-1) / k(1:ng-1)
end do
if (ncomp .eq. 1) then
! In the one component case the OZ solution is trivial.
ck(:, 1) = hk(:, 1) / (1.0_dp + rho(1)*hk(:, 1))
else ! Multicomponent OZ solution
! As above set up a unit matrix, and the diagonal R matrix
rhomat = 0.0_dp
unita = 0.0_dp
do i = 1, ncomp
rhomat(i,i) = rho(i)
unita(i,i) = 1.0_dp
end do
! Do the matrix calculations for each wavevector k.
do ik = 1, ng-1
! Convert the reciprocal space functions into matrices.
do j = 1, ncomp
do i = 1, j
ij = i + j*(j-1)/2
hmat(i, j) = hk(ik, ij)
if (i.lt.j) hmat(j, i) = hmat(i, j)
end do
end do
! Construct A = I + H . R, and B = H
a = unita + matmul(hmat, rhomat)
b = hmat
! Solve A.X = B so that X = (I + H.R)^(-1) . H.
call axeqb_reduce(a, ncomp, b, ncomp, perm, irc)
if (irc.gt.0) then
return_code = AXEQB_ERROR
error_msg = 'axeqb encountered singular problem in oz_solve2'
if (.not.suppress_msgs .and. .not.silent) then
print *, '** error: ', error_msg
print *, '** further messages of this kind will be suppressed'
suppress_msgs = .true.
end if
return
end if
! Now compute C = (I + H.R)^(-1) . H + beta UL
! (map back to functions, and unravel the pivoting)
! The + beta UL is done afterwards.
do j = 1, ncomp
do i = 1, j
ij = i + j*(j-1)/2
ck(ik, ij) = b(perm(i), j)
end do
end do
end do
end if
! Do the Fourier back transforms.
do i = 1, nfnc
fftwx(1:ng-1) = k(1:ng-1) * ck(1:ng-1, i)
call dfftw_execute(plan)
c(1:ng-1, i, 1) = (deltak / twopi**2) * fftwy(1:ng-1) / r(1:ng-1)
end do
! Add + beta UL (since we have the exact real and reciprocal space
! functions we don't need to pass this through the DDFT).
! Strictly speaking this is a flourish that isn't necessary for
! the structural features (pair functions and structure factors)
! since it cancels out again. However several of the
! thermodynamic calculations assume this offset is present.
ck = ck + ulongk
c(:,:,1) = c(:,:,1) + ulong
! Recover the indirect correlation functions. This means the
! results can be used in the structure and thermodynamics routines
! as though they had come from the HNC/MSA/RPA solution.
ek = hk - ck
e(:, :, 1) = h0(:, :) - c(:, :, 1)
end subroutine oz_solve2
! This routine implements the HNC closure expressed as c = exp(-beta v
! + e) - e - 1 where e = h - c is the indirect correlation function, c
! is the direct correlation function from the Ornstein-Zernicke
! relation, h = g - 1, and g is the pair distribution function. One
! can show this is equivalent to g = exp(-v + h - c) in Hansen +
! McDonald. As above, the routine actually works with c' = c + Ulong
! and e' = e - Ulong where Ulong is the long-range part of the
! potential for which the Fourier transform is simple. This means
! that 'v' in the above expression is the short-range part of the
! potential only. As we iterate we move forward in the history
! trajectory, cyclically.
subroutine hnc_picard
implicit none
integer :: i1, i0, i
istep = istep + 1
i1 = mod(istep-1, nps) + 1
i0 = i1 - 1; if (i0.eq.0) i0 = nps
do i = 1, nfnc
c(:,i,i1) = alpha * ( expnegus(:,i) * exp(e(:,i,i0)) &
& - e(:,i,i0) - 1.0_dp ) &
& + (1.0_dp - alpha) * c(:,i,i0)
end do
end subroutine hnc_picard
! The next routine implements the Ng method [K-C Ng, J. Chem. Phys.
! v61, 2680 (1974)] as an accelerated solver for the HNC closure. As
! we iterate we move forward in the history trajectory, cyclically.
subroutine hnc_ng
implicit none
integer :: i1, i0, j, j1, j2, p, nd, icp
real(kind=dp) :: dc(ng-1,nfnc,nps-1), de(ng-1,nfnc,nps-1), &
& a(nps-1,nps-1), x(nps-1), y(nps-1), yy, aux
integer :: ipiv(nps-1), info
integer, parameter :: lwork = 100