/
MultipleBranches.jl
687 lines (550 loc) · 26.6 KB
/
MultipleBranches.jl
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
# MultipleBranches.jl
# Extending the Feynman theory to multiple phonon branches
"""
frohlichPartial((f, ϵ_mode); ϵ_o,ϵ_s,meff)
Calculate a (partial) dielectric electron-phonon coupling element.
f - frequency of mode in THz
ϵ_mode - this mode's contribution to dielectric
ϵ_o - optical dielectric
ϵ_s - total static dielectric contribution
"""
function frohlichPartial((f, ϵ_mode); ϵ_o,ϵ_s,meff)
ω=f*1E12*2π
α=1/(4*π*ϵ_0) * ϵ_mode/(ϵ_o*ϵ_s) * (q^2/ħ) * sqrt(meff*me/(2*ω*ħ))
return α
end
# deprecated signature wrapped via multiple dispatch
frohlichPartial(ϵ_o,ϵ_s,ϵ_mode,f,meff) = frohlichPartial((f, ϵ_mode), ϵ_o=ϵ_o,ϵ_s=ϵ_s,meff=meff)
rows(M::Matrix) = map(x->reshape(getindex(M, x, :), :, size(M)[2]), 1:size(M)[1])
"""
IRtoDielectric(IRmodes,volume)
From absolute value of IR activities of phonon modes, generate a per-mode
contribution to the low-frequency dielectric constant.
IRmodes are tuples f,S with Frequency in THz; InfraRed activity in e^2 amu^-1
"""
function IRtoDielectric(IRmodes,volume)
ϵ=0.0 #* q^2 * THz^-2 * amu^-1 * m^-3
for r in rows(IRmodes)
f,S=r # frequency in THz; activity in e^2 amu^-1
f=f * 1E12 #* THz
ω=2π * f
S=S * q^2 / amu
ϵ_mode = S / ω^2 / volume
ϵ_mode /= 3 # take isotropic average = divide by 3
ϵ+=ϵ_mode
println("Mode f= $f S= $S ϵ_mode = $(ϵ_mode/ϵ_0)")
end
println("Raw ionic dielectric contribution: $ϵ absolute $(ϵ/ϵ_0) relative")
return ϵ/ϵ_0
end
"""
IRtoalpha(IR,volume)
Calculates contribution to dielectric constant for each polar phonon mode, and
thereby the Frohlich alpha contribution for this mode.
"""
function IRtoalpha(IR; volume, ϵ_o,ϵ_s,meff)
ϵ=0.0 #* q^2 * THz^-2 * amu^-1 * m^-3
α_sum=0.0
for r in rows(IR)
f,S=r # frequency in THz; activity in e^2 amu^-1
f=f * 1E12 #* THz
ω=2π * f
S=S * q^2 / amu
ϵ_mode = S / ω^2 / volume
ϵ_mode /= 3 # take isotropic average = divide by 3
ϵ_mode /= ϵ_0 # reduced dielectric constant
ϵ+=ϵ_mode
#println("Mode f= $f S= $S ϵ_mode = $(upreferred(ϵ_mode/u"ϵ0"))")
α=frohlichPartial(ϵ_o,ϵ_s,ϵ_mode,f/1E12,meff)
α_sum+=α
if (α>0.1)
println("Notable Mode f= $f α_partial=$α")
end
end
println("Sum alpha: $(α_sum)")
return α_sum
end
"""
DieletricFromIRmode(IRmode)
Calculate dielectric from an individual mode.
IRmode is a tuple f,S with Frequency in THz; InfraRed activity in e^2 amu^-1
"""
function DielectricFromIRmode(IRmode; volume)
f,S=IRmode # Assumes Frequency in THz; InfraRed activity in e^2 amu^-1
f*=1E12 # convert to Hz from THz
ω=2π * f
S=S * q^2 / amu
ϵ_mode = S / ω^2 / volume
ϵ_mode /= 3 # take isotropic average = divide by 3
ϵ_mode /= ϵ_0 # reduced dielectric
return ϵ_mode
end
function Hellwarth1999mobilityRHS( (α, (v,w) ,f) , effectivemass, T)
mb=effectivemass*MassElectron
ω=f*1E12*2π
βred=ħ*ω/(kB*T)
R=(v^2-w^2)/(w^2*v) # inline, page 300 just after Eqn (2)
b=R*βred/sinh(βred*v/2) # Feynman1962 version; page 1010, Eqn (47b)
a=sqrt( (βred/2)^2 + R*βred*coth(βred*v/2))
k(u,a,b,v) = (u^2+a^2-b*cos(v*u))^(-3/2)*cos(u) # integrand in (2)
K=quadgk(u->k(u,a,b,v),0,Inf)[1] # numerical quadrature integration of (2)
#Right-hand-side of Eqn 1 in Hellwarth 1999 // Eqn (4) in Baggio1997
RHS=α/(3*sqrt(π)) * βred^(5/2) / sinh(βred/2) * (v^3/w^3) * K
μ=RHS^-1 * (q)/(ω*mb)
return 1/μ
end
"""
----------------------------------------------------------------------
Multiple Branch Frohlich Alpha
----------------------------------------------------------------------
Partial dielectric electron-phonon coupling parameter, decomposed into ionic dielectric
contributions from each phonon mode of the matieral.
"""
"""
ϵ_ionic_mode(phonon_mode_freq::Float64, ir_activity::Float64, volume::Float64)
Calculate the ionic contribution to the dielectric function for a given phonon mode.
phonon_mode_freq is the frequency of the mode in THz.
- ir_activity is the infra-red activity of the mode in e^2 amu^-1.
- volume is the volume of the unit cell of the material in m^3.
"""
function ϵ_ionic_mode(phonon_mode_freq, ir_activity, volume) # single ionic mode
# Angular phonon frequency for the phonon mode (rad Hz)
ω_j = 2π * phonon_mode_freq * 1e12
# Dielectric contribution from a single ionic phonon mode
ϵ_mode = eV^2 * ir_activity / (3 * volume * ω_j^2 * amu)
# Normalise ionic dielectric contribution with 1 / (4π ϵ_0) (NB: the 4π has been pre-cancelled)
return ϵ_mode / ϵ_0
end
"""
ϵ_total(freqs_and_ir_activity::Matrix{Float64}, volume::Float64)
Calculate the total ionic contribution to the dielectric function from all phonon modes.
- freqs_and_ir_activity is a matrix containexing the phonon mode frequencies (in THz)
in the first column and the infra-red activities (in e^2 amu^-1) in the second column.
- volume is the volume of the unit cell of the material in m^3.
"""
function ϵ_total(freqs_and_ir_activity, volume) # total ionic contribution to dielectric
# Extract phonon frequencies (THz)
phonon_freqs = freqs_and_ir_activity[:, 1]
# Extra infra-red activities (e^2 amu^-1)
ir_activity = freqs_and_ir_activity[:, 2]
# Sum over all ionic contribution from each phonon mode
total_ionic = 0.0
for t in 1:length(phonon_freqs)
total_ionic += ϵ_ionic_mode(phonon_freqs[t], ir_activity[t], volume)
end
return total_ionic
end
"""
effective_freqs(freqs_and_ir_activity::Matrix{Float64}, num_var_params::Integer)
Generates a matrix of effective phonon modes with frequencies and infra-red activities
derived from a larger matrix using the Principal Component Analysis (PCA) method.
- freqs_and_ir_activity: is a matrix containing the phonon mode frequencies (in THz) in
the first column and the infra-red activities (in e^2 amu^-1) in the second column.
- num_var_params: is the number of effective modes required (which needs to be less
than the number of modes in freqs_and_ir_activity)
*** POSSIBLY REDUNDANT ***
"""
function effective_freqs(freqs_and_ir_activity, num_var_params) # PCA Algorithm
# Check that the number of effective modes is less than the number of actual phonon modes.
if num_var_params >= size(freqs_and_ir_activity)[1]
println("The number of effective phonon modes has to be less than the total number of phonon modes.")
else
# Centralise data by subtracting the columnwise mean
standardized_matrix = freqs_and_ir_activity' .- mean(freqs_and_ir_activity', dims = 2)
# Calculate the covariance matrix S' * S. Matrix size is (n - 1) x (n - 1) for number of params (here n = 2)
covariance_matrix = standardized_matrix' * standardized_matrix
# Extract eigenvectors of the covariance matrix
eigenvectors = eigvecs(covariance_matrix)
# Project the original data along the covariance matrix eigenvectors and undo the centralisation
reduced_matrix = standardized_matrix[:, 1:num_var_params] * eigenvectors[1:num_var_params, 1:num_var_params] *
eigenvectors[1:num_var_params, 1:num_var_params]' .+ mean(freqs_and_ir_activity', dims = 2)
# Resultant matrix is positive definite and transposed.
return abs.(reduced_matrix')
end
end
"""
multi_frohlichalpha(ϵ_optic::Float64, ϵ_ionic::Float64, ϵ_total::Float64, phonon_mode_freq::Float64, m_eff::Float64)
Calculates the partial dielectric electron-phonon coupling parameter for a given
longitudinal optical phonon mode. This decomposes the original Frohlich alpha coupling
parameter (defined for a single phonon branch) into contributions from multiple phonon
branches.
- ϵ_optic is the optical dielectric constant of the material.
- ϵ_ionic is the ionic dielectric contribution from the phonon mode.
- ϵ_total is the total ionic dielectric contribution from all phonon modes of the material.
- phonon_mode_freq is the frequency of the phonon mode (THz).
- m_eff is the band mass of the electron (in units of electron mass m_e)
"""
function multi_frohlichalpha(ϵ_optic, ϵ_ionic, ϵ_total, phonon_mode_freq, m_eff)
# The Rydberg energy unit
Ry = eV^4 * me / (2 * ħ^2)
# Angular phonon frequency for the phonon mode (rad Hz).
ω = 2π * 1e12 * phonon_mode_freq
# The static dielectric constant. Calculated here instead of inputted so that ionic modes are properly normalised.
ϵ_static = ϵ_total + ϵ_optic
# The contribution to the electron-phonon parameter from the currrent phonon mode. 1 / (4π ϵ_0) is the dielectric normalisation.
α_j = (m_eff * Ry / (ħ * ω))^(1 / 2) * ϵ_ionic / (4π * ϵ_0) / (ϵ_optic * ϵ_static)
return α_j
end
"""
----------------------------------------------------------------------
Multiple Branch Polaron Free Energy
----------------------------------------------------------------------
Calculate the polaron free energy, generalised from Osaka's expression to the case where
multiple phonon modes are present in the material.
"""
"""
κ_i(i::Integer, v::Array{Float64}(undef, 1), w::Array{Float64}(undef, 1))
Calculates the spring-constant coupling the electron to the `ith' fictitious mass that
approximates the exact electron-phonon interaction with a harmonic coupling to a massive
fictitious particle. NB: Not to be confused with the number of physical phonon branches;
many phonon branches could be approximated with one or two etc. fictitious masses for
example. The number of fictitious mass does not necessarily need to match the number of
phonon branches.
- i enumerates the current fictitious mass.
- v is an one-dimensional array of the v variational parameters.
- w is an one-dimensional array of the w variational parameters.
"""
function κ_i(i, v, w)
κ = v[i]^2 - w[i]^2
if length(v) > 1
for j in 1:length(v)
if j != i
κ *= (v[j]^2 - w[i]^2) / (w[j]^2 - w[i]^2)
end
end
end
return κ
end
"""
h_i(i::Integer, v::Array{Float64}(undef, 1), w::Array{Float64}(undef, 1))
Calculates the normal-mode (the eigenmodes) frequency of the coupling between the
electron and the `ith' fictitious mass that approximates the exact electron-phonon
interaction with a harmonic coupling to a massive fictitious particle. NB: Not to be
confused with the number of physical phonon branches; many phonon branches could be
approximated with one or two etc. fictitious masses for example. The number of
fictitious mass does not necessarily need to match the number of phonon branches.
- i enumerates the current fictitious mass.
- v is an one-dimensional array of the v variational parameters.
- w is an one-dimensional array of the w variational parameters.
"""
function h_i(i, v, w)
h = v[i]^2 - w[i]^2
if length(v) > 1
for j in 1:length(v)
if j != i
h *= (w[j]^2 - v[i]^2) / (v[j]^2 - v[i]^2)
end
end
end
return h
end
"""
C_ij(i::Integer, j::Integer, v::Array{Float64}(undef, 1), w::Array{Float64}(undef, 1))
Calculates the element to the coupling matrix C_ij (a generalisation of Feynman's `C`
coupling variational parameter) between the electron and the `ith' and `jth' fictitious
masses that approximates the exact electron-phonon interaction with a harmonic coupling
to a massive fictitious particle. NB: Not to be confused with the number of physical
phonon branches; many phonon branches could be approximated with one or two etc.
fictitious masses for example. The number of fictitious mass does not necessarily need
to match the number of phonon branches.
- i, j enumerate the current fictitious masses under focus (also the index of the element in the coupling matrix C)
- v is an one-dimensional array of the v variational parameters.
- w is an one-dimensional array of the w variational parameters.
"""
function C_ij(i, j, v, w)
C = w[i] * κ_i(i, v, w) * h_i(j, v, w) / (4 * (v[j]^2 - w[i]^2))
return C
end
"""
D_j(τ::Float64, β::Float64, v::Array{Float64}(undef, 1), w::Array{Float64}(undef, 1))
Calculates the recoil function that approximates the exact influence (recoil effects) of
the phonon bath on the electron with the influence of the harmonicly coupled fictitious
masses on the electron. It appears in the exponent of the intermediate scattering
function.
- τ is the imaginary time variable.
- β is the reduced thermodynamic temperature ħ ω_j / (kB T) associated with the `jth` phonon mode.
- v is an one-dimensional array of the v variational parameters.
- w is an one-dimensional array of the w variational parameters.
"""
function D_j(τ, β, v, w)
D = τ * (1 - τ / β)
for i in 1:length(v)
D += (h_i(i, v, w) / v[i]^2) * ((1 + exp(-v[i] * β) - exp(-v[i] * τ) - exp(v[i] * (τ - β))) / (v[i] * (1 - exp(-v[i] * β))) - τ * (1 - τ / β))
end
return D
end
"""
D_j(τ::Float64, v::Array{Float64}(undef, 1), w::Array{Float64}(undef, 1))
Calculates the recoil function that approximates the exact influence (recoil effects) of
the phonon bath on the electron with the influence of the harmonicly coupled fictitious
masses on the electron. It appears in the exponent of the intermediate scattering
function. This function works at zero temperature.
- τ is the imaginary time variable.
- v is an one-dimensional array of the v variational parameters.
- w is an one-dimensional array of the w variational parameters.
"""
function D_j(τ, v, w)
D = τ
for i in 1:length(v)
D += (h_i(i, v, w) / v[i]^2) * ((1 - exp(-v[i] * τ)) / v[i] - τ)
end
return D
end
"""
B_j(α::Float64, β::Float64, v::Array{Float64}(undef, 1), w::Array{Float64}(undef, 1))
Generalisation of the B function from Equation 62c in Hellwarth and Biaggio[1]. This is
the expected value of the exact action <S_j> taken w.r.t trial action, given for the
`jth` phonon mode.
- α is the partial dielectric electron-phonon coupling parameter for the `jth` phonon
mode.
- β is the reduced thermodynamic temperature ħ ω_j / (kB T) associated with the
`jth` phonon mode.
- v is an one-dimensional array of the v variational parameters.
- w is an one-dimensional array of the w variational parameters.
[1] Hellwarth, R. W., Biaggio, I. (1999). Mobility of an electron in a multimode polar
lattice. Physical Review B, 60(1), 299–307. https://doi.org/10.1103/PhysRevB.60.299
"""
function B_j(α, β, v, w; rtol = 1e-3)
B_integrand(τ) = (exp(-τ) + exp(τ - β)) / sqrt(abs(D_j(τ, β, v, w)))
B = α / (√π * (1 - exp(-β))) * quadgk(τ -> B_integrand(τ), 0.0, β / 2, rtol = rtol)[1]
return B
end
"""
B_j(α::Float64, v::Array{Float64}(undef, 1), w::Array{Float64}(undef, 1))
Generalisation of the B function from Equation 62c in Biaggio and Hellwarth [1] to zero
temperature. This is the expected value of the exact action <S_j> taken w.r.t trial
action, given for the `jth` phonon mode.
- α is the partial dielectric electron-phonon coupling parameter for the `jth` phonon mode.
- v is an one-dimensional array of the v variational parameters.
- w is an one-dimensional array of the w variational parameters.
[1] Hellwarth, R. W., Biaggio, I. (1999). Mobility of an electron in a multimode polar
lattice. Physical Review B, 60(1), 299–307. https://doi.org/10.1103/PhysRevB.60.299
"""
function B_j(α, v, w; rtol = 1e-3)
B_integrand(τ) = exp(-abs(τ)) / sqrt(abs(D_j(abs(τ), v, w)))
B = α / √π * quadgk(τ -> B_integrand(τ), 0.0, Inf64, rtol = rtol)[1]
return B
end
"""
C_j(β::Float64, v::Array{Float64}(undef, 1), w::Array{Float64}(undef, 1), n::Float64)
Generalisation of the C function from Equation 62e in Hellwarth and Biaggio[1]. This is
the expected value of the trial action <S_0> taken w.r.t trial action.
- β is the reduced thermodynamic temperature ħ ω_j / (kB T) associated with the `jth`
phonon mode.
- v is an one-dimensional array of the v variational parameters.
- w is an one-dimensional array of the w variational parameters.
- n is the total number of phonon modes.
[1] Hellwarth, R. W., Biaggio, I. (1999). Mobility of an electron in a multimode polar
lattice. Physical Review B, 60(1), 299–307. https://doi.org/10.1103/PhysRevB.60.299
"""
function C_j(β, v, w, n)
s = 0.0
# Sum over the contributions from each fictitious mass.
for i in 1:length(v)
for j in 1:length(v)
s += C_ij(i, j, v, w) / (v[j] * w[i]) * (coth(β * v[j] / 2) - 2 / (β * v[j]))
end
end
# Divide by the number of phonon modes to give an average contribution per phonon mode.
return 3 * s / n
end
"""
C_j(v::Array{Float64}(undef, 1), w::Array{Float64}(undef, 1), n::Float64)
Generalisation of the C function from Equation 62e in Biaggio and Hellwarth [] but to
zero temperaure. This is the expected value of the trial action <S_0> taken w.r.t trial
action.
- v is an one-dimensional array of the v variational parameters.
- w is an one-dimensional array of the w variational parameters.
- n is the total number of phonon modes.
[1] Hellwarth, R. W., Biaggio, I. (1999). Mobility of an electron in a multimode polar
lattice. Physical Review B, 60(1), 299–307. https://doi.org/10.1103/PhysRevB.60.299
"""
function C_j(v, w, n)
s = 0.0
# Sum over the contributions from each fictitious mass.
for i in 1:length(v)
for j in 1:length(v)
s += C_ij(i, j, v, w) / (v[j] * w[i])
end
end
# Divide by the number of phonon modes to give an average contribution per phonon mode.
return 3 * s / n
end
"""
A_j(β::Float64, v::Array{Float64}(undef, 1), w::Array{Float64}(undef, 1), n::Float64)
Generalisation of the A function from Equation 62b in Hellwarth and Biaggio[1]. This is
the Helmholtz free energy of the trial model.
- β is the reduced thermodynamic temperature ħ ω_j / (kB T) associated with the `jth` phonon mode.
- v is an one-dimensional array of the v variational parameters.
- w is an one-dimensional array of the w variational parameters.
- n is the total number of phonon modes.
[1] Hellwarth, R. W., Biaggio, I. (1999). Mobility of an electron in a multimode polar
lattice. Physical Review B, 60(1), 299–307. https://doi.org/10.1103/PhysRevB.60.299
"""
function A_j(β, v, w, n)
s = -log(2π * β) / 2
# Sum over the contributions from each fictitious mass.
for i in 1:length(v)
if v[i] != w[i]
s += v[i] == w[i] ? 0.0 : log(v[i] / w[i]) - log(sinh(v[i] * β / 2) / sinh(w[i] * β / 2))
end
end
# Divide by the number of phonon modes to give an average contribution per phonon mode.
3 / β * s / n
end
"""
A_j(v::Array{Float64}(undef, 1), w::Array{Float64}(undef, 1), n::Float64)
Generalisation of the A function from Equation 62b in Hellwarth Biaggio[1] but to
zero temperature. This is the ground-state energy of the trial model.
- v is an one-dimensional array of the v variational parameters.
- w is an one-dimensional array of the w variational parameters.
- n is the total number of phonon modes.
[1] Hellwarth, R. W., Biaggio, I. (1999). Mobility of an electron in a multimode polar
lattice. Physical Review B, 60(1), 299–307. https://doi.org/10.1103/PhysRevB.60.299
"""
function A_j(v, w, n)
s = sum(v .- w)
return -3 * s / (2 * n)
end
"""
multi_F(v, w, α, β; ω = 1.0)
Calculates the Helmholtz free energy of the polaron for a material with multiple phonon
branches. This generalises Osaka's free energy expression (below Equation (22) in [1]).
- v is an one-dimensional array of the v variational parameters.
- w is an one-dimensional array of the w variational parameters.
- α is the Frohlich coupling constant.
- β is the thermodynamic temperature.
[1] Osaka, Y. (1959). Polaron State at a Finite Temperature. Progress of Theoretical
Physics, 22(3), 437–446. https://doi.org/10.1143/ptp.22.437
"""
function multi_F(v, w, α, β; ω = 1.0, rtol = 1e-3, T = nothing, verbose = false)
num_modes = length(ω)
# Add contribution to the total free energy from the phonon mode.
F = sum(-(B_j(α[j], β[j], v, w; rtol = rtol) + C_j(β[j], v, w, num_modes) + A_j(β[j], v, w, num_modes)) * ω[j] for j in 1:num_modes)
# Print the free energy.
if verbose
println("\e[2K", "Process: $(count) / $processes ($(round.(count / processes * 100, digits = 1)) %) | T = $(round.(T, digits = 3)) | F = $(round.(F, digits = 3))")
print("\033[F")
global count += 1
end
# Free energy in units of meV
return F
end
"""
multi_F(v, w, α; ω)
Calculates the Helmholtz free energy of the polaron for a material with multiple phonon
branches. This generalises Osaka's free energy expression (below Equation (22) in [1]).
- α is the Frohlich alpha parameter.
- ω is a vector containing the phonon mode frequencies (in THz).
- v and w determines if the function should start with a random initial set of
variational parameters (= 0.0) or a given set of variational parameter values.
[1] Osaka, Y. (1959). Polaron State at a Finite Temperature. Progress of Theoretical
Physics, 22(3), 437–446. https://doi.org/10.1143/ptp.22.437
"""
function multi_F(v, w, α; ω = 1.0, rtol = 1e-3, verbose = false)
num_modes = length(ω)
# Add contribution to the total free energy from the phonon mode.
F = sum(-(B_j(α[j], v, w; rtol = rtol) + C_j(v, w, num_modes) + A_j(v, w, num_modes)) * ω[j] for j in 1:num_modes)
# Print the free energy.
if verbose
println("\e[2K", "Process: $(count) / $processes ($(round.(count / processes * 100, digits = 1)) %) | T = 0.0 | F = $(round.(F, digits = 3))")
print("\033[F")
global count += 1
end
# Free energy in units of meV
return F
end
"""
variation(α::Vector{Real}, β::Vector{Real}; v::Real, w::Real, ω::Vector{Real}, N::Integer)
Minimises the multiple phonon mode free energy function for a set of v_p and w_p
variational parameters. The variational parameters follow the inequality: v_1 > w_1
> v_2 > w_2 > ... > v_N > w_N.
- β is the thermodynamic temperature.
- α is the Frohlich alpha parameter.
- ω is a vecotr containing the phonon mode frequencies (in THz).
- v and w determines if the function should start with a random initial set of variational parameters (= 0.0) or a given set of variational parameter values.
- N specifies the number of variational parameter pairs, v_p and w_p, to use in minimising the free energy.
"""
function var_params(α, β; v = 0.0, w = 0.0, ω = 1.0, N = 1, rtol = 1e-3, show_trace = false, T = nothing, verbose = false) # N number of v and w params
if N != length(v) != length(w)
return error("The number of variational parameters v & w must be equal to N.")
end
# Use a random set of N initial v and w values.
if v == 0.0 || w == 0.0
# Intial guess for v and w parameters.
initial = sort(rand(2 * N), rev=true) .* 4.0 .+ 1.0 # initial guess around 4 and ≥ 1.
else
Δv = v .- w
initial = vcat(Δv .+ rtol, w)
end
# Limits of the optimisation.
lower = fill(0.0, 2 * N)
upper = fill(1000.0, 2 * N)
# The multiple phonon mode free energy function to minimise.
f(x) = multi_F([x[2 * n - 1] for n in 1:N] .+ [x[2 * n] for n in 1:N], [x[2 * n] for n in 1:N], α, β; ω = ω, rtol = rtol)
# Use Optim to optimise the free energy function w.r.t the set of v and w parameters.
solution = Optim.optimize(
Optim.OnceDifferentiable(f, initial; autodiff = :forward),
lower,
upper,
initial,
Fminbox(BFGS()),
Optim.Options(show_trace = show_trace), # Set time limit for asymptotic convergence if needed.
)
# Extract the v and w parameters that minimised the free energy.
var_params = Optim.minimizer(solution)
if Optim.converged(solution) == false
@warn "Failed to converge feynmanvw solution."
end
# Separate the v and w parameters into one-dimensional arrays (vectors).
Δv = [var_params[2 * n - 1] for n in 1:N]
w = [var_params[2 * n] for n in 1:N]
# Print the variational parameters that minimised the free energy.
if verbose
println("\e[2K", "Process: $(count) / $processes ($(round.(count / processes * 100, digits = 1)) %) | T = $(round.(T, digits = 3)) | v = $(round.(Δv .+ w, digits = 3)) | w = $(round.(w, digits = 3))")
print("\033[F")
global count += 1
end
# Return the variational parameters that minimised the free energy.
return hcat(Δv .+ w, w)
end
function var_params(α; v = 0.0, w = 0.0, ω = 1.0, N = 1, rtol = 1e-3, show_trace = false, verbose = false) # N number of v and w params
if N != length(v) != length(w)
return error("The number of variational parameters v & w must be equal to N.")
end
# Use a random set of N initial v and w values.
if v == 0.0 || w == 0.0
# Intial guess for v and w parameters.
initial = sort(rand(2 * N), rev=true) .* 4.0 .+ 1.0 # initial guess around 4 and ≥ 1.
else
Δv = v .- w
initial = vcat(Δv .+ rtol, w)
end
# Limits of the optimisation.
lower = fill(0.0, 2 * N)
upper = fill(1000, 2 * N)
# The multiple phonon mode free energy function to minimise.
f(x) = multi_F([x[2 * n - 1] for n in 1:N] .+ [x[2 * n] for n in 1:N], [x[2 * n] for n in 1:N], α; ω = ω, rtol = rtol)
# Use Optim to optimise the free energy function w.r.t the set of v and w parameters.
solution = Optim.optimize(
Optim.OnceDifferentiable(f, initial; autodiff = :forward),
lower,
upper,
initial,
Fminbox(BFGS()),
Optim.Options(show_trace = show_trace), # Set time limit for asymptotic convergence if needed.
)
# Extract the v and w parameters that minimised the free energy.
var_params = Optim.minimizer(solution)
# Separate the v and w parameters into one-dimensional arrays (vectors).
Δv = [var_params[2 * n - 1] for n in 1:N]
w = [var_params[2 * n] for n in 1:N]
# Print the variational parameters that minimised the free energy.
if verbose
println("\e[2K", "Process: $(count) / $processes ($(round.(count / processes * 100, digits = 1)) %) | T = 0.0 | v = $(round.(Δv .+ w, digits = 3)) | w = $(round.(w, digits = 3))")
print("\033[F")
global count += 1
end
# Return the variational parameters that minimised the free energy.
return hcat(Δv .+ w, w)
end