-
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathpartition_functions.f90
566 lines (515 loc) · 26.1 KB
/
partition_functions.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
!=========================================================================================
! Peacemaker -- A Quantum Cluster Equilibrium Code.
!
! Copyright 2004-2006 Barbara Kirchner, University of Bonn
! Copyright 2007-2012 Barbara Kirchner, University of Leipzig
! Copyright 2013-2018 Barbara Kirchner, University of Bonn
!
! This file is part of Peacemaker.
!
! Peacemaker 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.
!
! Peacemaker 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 Peacemaker. If not, see <http://www.gnu.org/licenses/>
!=========================================================================================
! This module implements the partition functions used by Peacemaker.
module partition_functions
use kinds
use cluster
use constants
use shared_data
implicit none
private
!=====================================================================================
! Public entities
public :: calculate_lnq
public :: update_lnq
public :: calculate_dlnq
public :: calculate_ddlnq
public :: pf_t
!=====================================================================================
! Data type representing a cluster partition function.
type :: pf_t
real(dp) :: qtrans, qvib, qrot, qelec, qint, qtot
end type pf_t
!=====================================================================================
contains
!=================================================================================
! Calculates all cluster partition functions.
subroutine calculate_lnq(lnq, amf, bxv, temp, vol)
type(pf_t), dimension(size(clusterset)), intent(out) :: lnq
real(dp), intent(in) :: amf
real(dp), intent(in) :: bxv
real(dp), intent(in) :: vol
real(dp), intent(in) :: temp
call calculate_lnqtrans(lnq(:)%qtrans, bxv, temp, vol)
call calculate_lnqvib(lnq(:)%qvib, temp)
call calculate_lnqrot(lnq(:)%qrot, temp)
call calculate_lnqelec(lnq(:)%qelec, temp)
call calculate_lnqint(lnq(:)%qint, amf, temp, vol)
lnq(:)%qtot = lnq(:)%qtrans + lnq(:)%qvib + lnq(:)%qrot + lnq(:)%qelec + &
lnq(:)%qint
end subroutine calculate_lnq
!=================================================================================
! Updates the cluster partition functions that depend on the volume.
subroutine update_lnq(lnq, amf, bxv, temp, vol)
type(pf_t), dimension(size(clusterset)), intent(inout) :: lnq
real(dp), intent(in) :: amf
real(dp), intent(in) :: bxv
real(dp), intent(in) :: temp
real(dp), intent(in) :: vol
call calculate_lnqtrans(lnq(:)%qtrans, bxv, temp, vol)
call calculate_lnqint(lnq(:)%qint, amf, temp, vol)
lnq(:)%qtot = lnq(:)%qtrans + lnq(:)%qvib + lnq(:)%qrot + lnq(:)%qelec + &
lnq(:)%qint
end subroutine update_lnq
!=================================================================================
! Calculates analytical derivatives of the cluster partition functions.
subroutine calculate_dlnq(dlnq, amf, bxv, amf_temp, bxv_temp, temp, vol)
type(pf_t), dimension(size(clusterset)), intent(out) :: dlnq
real(dp), intent(in) :: amf
real(dp), intent(in) :: bxv
real(dp), intent(in) :: amf_temp
real(dp), intent(in) :: bxv_temp
real(dp), intent(in) :: temp
real(dp), intent(in) :: vol
call calculate_dlnqtrans(dlnq(:)%qtrans, bxv, bxv_temp, temp, vol)
call calculate_dlnqvib(dlnq(:)%qvib, temp)
call calculate_dlnqrot(dlnq(:)%qrot, temp)
call calculate_dlnqelec(dlnq(:)%qelec, temp)
call calculate_dlnqint(dlnq(:)%qint, temp, amf, amf_temp, vol)
dlnq(:)%qtot = dlnq(:)%qtrans + dlnq(:)%qvib + dlnq(:)%qrot + &
dlnq(:)%qelec + dlnq(:)%qint
end subroutine calculate_dlnq
!=================================================================================
! Calculates analytical curvatures of the cluster partition functions.
subroutine calculate_ddlnq(ddlnq, amf, bxv, amf_temp, bxv_temp, temp, vol)
type(pf_t), dimension(size(clusterset)), intent(out) :: ddlnq
real(dp), intent(in) :: amf
real(dp), intent(in) :: bxv
real(dp), intent(in) :: amf_temp
real(dp), intent(in) :: bxv_temp
real(dp), intent(in) :: temp
real(dp), intent(in) :: vol
call calculate_ddlnqtrans(ddlnq(:)%qtrans, bxv, bxv_temp, temp, vol)
call calculate_ddlnqvib(ddlnq(:)%qvib, temp)
call calculate_ddlnqrot(ddlnq(:)%qrot, temp)
call calculate_ddlnqelec(ddlnq(:)%qelec, temp)
call calculate_ddlnqint(ddlnq(:)%qint, temp, amf, amf_temp, vol)
ddlnq(:)%qtot = ddlnq(:)%qtrans + ddlnq(:)%qvib + ddlnq(:)%qrot + &
ddlnq(:)%qelec + ddlnq(:)%qint
end subroutine calculate_ddlnq
!=================================================================================
! Calculates the translational cluster partition function.
! q_trans = (2*pi*m*kb*T/h^2)^1.5*V
subroutine calculate_lnqtrans(lnq, bxv, temp, vol)
real(dp), dimension(:), intent(out) :: lnq
real(dp), intent(in) :: vol
real(dp), intent(in) :: temp
real(dp), intent(in) :: bxv
real(dp):: lambda
real(dp):: mass
integer:: iclust
do iclust = 1, size(clusterset)
! Calculate the de Broglie wave length
mass = clusterset(iclust)%mass*amu
lambda = planck/sqrt(2.0_dp*pi*mass*kb*temp)
! Calculate the partition function
lnq(iclust) = log(vol-bxv*global_data%vexcl) - 3.0_dp*log(lambda)
end do
end subroutine calculate_lnqtrans
!=================================================================================
! Calculates the vibrational cluster partition function.
! edit: q_vib is the product, the sum occurs due to the logarithm
! q_vib = product over all nu: exp[-h*nu/(2*kb*T)]/{1-exp[-h*nu/(kb*T)]}
subroutine calculate_lnqvib(lnq, temp)
real(dp), dimension(:), intent(out) :: lnq
real(dp), intent(in) :: temp
integer:: iclust
integer:: ifreq, imoment
real(dp):: factor
real(dp):: t_vib
real(dp):: q_ho, q_fr
real(dp):: moi, emoi, Bav, t_rot, w
logical :: free_rotator
free_rotator = .false.
w = 1.0_dp
if (global_data%rotor_cutoff > 0.0_dp) free_rotator = .true.
! TODO: Check atoms.
factor = planck*100.0_dp*speed_of_light/kb
do iclust = 1, size(clusterset)
associate(f => clusterset(iclust)%frequencies, q => lnq(iclust), &
x => clusterset(iclust)%anharmonicity, c => clusterset(iclust))
q = 0.0_dp
q_ho = 0.0_dp
q_fr = 0.0_dp
if (c%atom) cycle
! Average intertia moment of cluster for free rotator approximation
Bav = sum(c%inertia)/real(size(c%inertia,1),dp)*amu*1.0e-20_dp
! Bav = 1.0_dp
! if (free_rotator) then
! do imoment = 1, size(c%inertia)
! Bav = Bav * c%inertia(imoment)*amu*1.0e-20_dp
! end do
! Bav = Bav**(1.0_dp/real(size(c%inertia), dp))
! end if
do ifreq = 1, size(f)
t_vib = factor*f(ifreq)
if (x > 0.0_dp) then
! Morse oscillator.
q = q - log(2.0_dp*sinh(0.5_dp*t_vib/temp)) + x*t_vib/temp * &
(0.25_dp + 0.5_dp/sinh(0.5_dp*t_vib/temp)**2)
else
! Harmonic oscillator.
! TODO: Use log1p and exp1m family of functions.
q_ho = -0.5_dp*t_vib/temp - log(1.0_dp - exp(-t_vib/temp))
if (free_rotator) then
! Free rotator approximation by Grimme.
moi = planck/(8 * pi**2 * 100.0_dp * speed_of_light * f(ifreq))
emoi = moi * Bav/(moi + Bav)
t_rot = hbar**2/(2.0_dp*emoi*kb)
q_fr = log(sqrt(pi*temp/t_rot)/real(c%sigma, dp))
w = 1.0_dp/(1.0_dp + (global_data%rotor_cutoff/f(ifreq))**4)
end if
q = q + w * q_ho + (1.0_dp - w) * q_fr
end if
end do
end associate
end do
end subroutine calculate_lnqvib
!=================================================================================
! Calculates the rotational cluster partition function.
! q_rot = (1/sigma)*(T/t_rot) --- linear
! q_rot = sqrt(pi)/sigma*(T^3/(t_rot1*t_rot2*t_rot3))^(1/2) --- nonlinear
! q_rot = 1 --- atom
subroutine calculate_lnqrot(lnq, temp)
real(dp), dimension(:), intent(out) :: lnq
real(dp), intent(in) :: temp
integer:: iclust
integer:: imoment
real(dp):: t_rot
do iclust = 1, size(clusterset)
associate(c => clusterset(iclust))
if (c%atom) then
lnq(iclust) = 0.0_dp
else
lnq(iclust) = 1.0_dp
do imoment = 1, size(c%inertia)
t_rot = hbar**2/(2.0_dp*c%inertia(imoment)*amu*1.0e-20_dp*kb)
lnq(iclust) = lnq(iclust) * temp/t_rot
end do
if (c%linear) then
! linear cluster:
! q_rot = T/(sigma*t_rot)
! Note: For a linear molecule, the inertia array contains
! two moments of inertia of equal magnitude, which were
! multiplied together above.
lnq(iclust) = sqrt(lnq(iclust))/real(c%sigma, dp)
else
! polyatomic, nonlinear cluster:
! q_rot = sqrt(pi)/sigma*(T^3/(t_rot1*t_rot2*t_rot3))
lnq(iclust) = sqrt(pi*lnq(iclust))/real(c%sigma, dp)
end if
! Let's not forget the logarithm.
lnq(iclust) = log(lnq(iclust))
end if
end associate
end do
end subroutine calculate_lnqrot
!=================================================================================
! Calculates the electronic cluster partition function.
! q_elec = exp(-dE/(kb*T))
subroutine calculate_lnqelec(lnq, temp)
real(dp), dimension(:), intent(out) :: lnq
real(dp), intent(in) :: temp
! The adiabatic interaction energy is in units of kJ/mol.
lnq(:) = (-1000.0_dp/avogadro)*clusterset(:)%energy/(kb*temp)
end subroutine calculate_lnqelec
!=================================================================================
! Calculates the mean field partition function.
subroutine calculate_lnqint(lnq, amf, temp, vol)
real(dp), dimension(:), intent(out) :: lnq
real(dp), intent(in) :: vol
real(dp), intent(in) :: amf
real(dp), intent(in) :: temp
integer:: iclust
real(dp):: emf
! do iclust = 1, size(clusterset)
! associate(c => clusterset(iclust))
! call calculate_emf(emf, c%composition, amf, vol)
! lnq(iclust) = -emf / (kb*temp)
! end associate
! end do
do iclust = 1, size(clusterset)
associate(c => clusterset(iclust))
emf = -amf*real(sum(c%composition), dp)*sum(global_data%ntot)/vol
lnq(iclust) = -emf / (kb*temp)
end associate
end do
end subroutine calculate_lnqint
!=================================================================================
! Calculates the mean field energy.
subroutine calculate_emf(emf, composition, amf, vol)
real(dp), intent(out) :: emf
integer, dimension(size(monomer)), intent(in) :: composition
real(dp), intent(in) :: amf, vol
integer :: i, j
real(dp), dimension(size(monomer)) :: amf_pure
! Choose old-style amf_pure, if not specified in input, i.e., if it is zero
if (all(global_data%amf_pure <= 0.0_dp)) then
amf_pure = amf
else
amf_pure = global_data%amf_pure
end if
! Calculate mean field energy
emf = 0.0_dp
do i = 1, size(monomer)
emf = emf + &
global_data%ntot(i)*real(composition(i), dp)*amf_pure(i)
do j = i+1, size(monomer)
emf = emf + &
global_data%ntot(i)*real(composition(j), dp)*amf + &
global_data%ntot(j)*real(composition(i), dp)*amf
end do
end do
emf = -1.0_dp/vol * emf
end subroutine calculate_emf
!=================================================================================
! Calculates the temperature derivative of the translational partition function.
subroutine calculate_dlnqtrans(dlnq, bxv, bxv_temp, temp, vol)
real(dp), dimension(:), intent(out) :: dlnq
real(dp), intent(in) :: bxv
real(dp), intent(in) :: bxv_temp
real(dp), intent(in) :: temp
real(dp), intent(in) :: vol
dlnq(:) = 1.5_dp/temp + &
(global_data%vexcl * bxv_temp)/(-vol+bxv*global_data%vexcl)
end subroutine calculate_dlnqtrans
!=================================================================================
! Calculates the temperature derivative of the vibrational partition function.
subroutine calculate_dlnqvib(dlnq, temp)
real(dp), dimension(:), intent(out) :: dlnq
real(dp), intent(in) :: temp
integer :: iclust
integer:: ifreq
real(dp):: factor
real(dp):: t_vib
real(dp):: fsinh
real(dp):: fcosh
real(dp):: q_ho, q_fr
real(dp):: w
! TODO: Check atoms.
factor = planck*100.0_dp*speed_of_light/kb
do iclust = 1, size(clusterset)
associate(f => clusterset(iclust)%frequencies, q => dlnq(iclust), &
x => clusterset(iclust)%anharmonicity)
q = 0.0_dp
q_ho = 0.0_dp
q_fr = 0.0_dp
do ifreq = 1, size(f)
t_vib = factor*f(ifreq)
if (x > 0.0_dp) then
! Morse oscillator.
fsinh = sinh(0.5_dp*t_vib/temp)
fcosh = cosh(0.5_dp*t_vib/temp)
q = q + 0.5_dp*t_vib*fcosh/(temp**2*fsinh) - &
x*t_vib/temp**2*(0.25_dp + 0.5_dp/fsinh**2) + &
0.5_dp*x*t_vib**2*fcosh/(temp**3*fsinh**3)
else
! Harmonic oscillator.
! TODO: Use log1p and exp1m family of functions.
q_ho = 0.5_dp*t_vib/temp**2 + &
t_vib/temp**2/(exp(t_vib/temp)-1.0_dp)
! Free rotator approximation by Grimme.
q_fr = 0.5_dp/temp
w = 1.0_dp/(1.0_dp + (global_data%rotor_cutoff/f(ifreq))**4)
q = q + w * q_ho + (1.0_dp - w) * q_fr
end if
end do
end associate
end do
end subroutine calculate_dlnqvib
!=================================================================================
! Calculates the temperature derivative of the rotational partition function.
subroutine calculate_dlnqrot(dlnq, temp)
real(dp), dimension(:), intent(out) :: dlnq
real(dp), intent(in) :: temp
integer :: iclust
do iclust = 1, size(clusterset)
associate(c => clusterset(iclust), q => dlnq(iclust))
if (c%atom) then
q = 0.0_dp
else if (c%linear) then
! linear cluster:
! q_rot = T/(sigma*t_rot) -> d/dt(lnq)=1/T
q = 1.0_dp/temp
else
! polyatomic, nonlinear cluster: d/dT(lnq)=3/(2T)
q = 1.5_dp/temp
end if
end associate
end do
end subroutine calculate_dlnqrot
!=================================================================================
! Calculates the temperature derivative of the electronic partition function.
subroutine calculate_dlnqelec(dlnq, temp)
real(dp), dimension(:), intent(out) :: dlnq
real(dp), intent(in) :: temp
! The adiabatic interaction energy is in units of kJ/mol.
dlnq(:) = (1000.0_dp/avogadro)*clusterset(:)%energy/(kb*temp**2)
end subroutine calculate_dlnqelec
!=================================================================================
! Calculates the temperature derivative of the mean field partition function.
subroutine calculate_dlnqint(dlnq, temp, amf, amf_temp, vol)
real(dp), dimension(:), intent(out) :: dlnq
real(dp), intent(in) :: temp
real(dp), intent(in) :: amf
real(dp), intent(in) :: amf_temp
real(dp), intent(in) :: vol
integer :: iclust
real(dp):: emf
! do iclust = 1, size(clusterset)
! associate(c => clusterset(iclust))
! call calculate_emf(emf, c%composition, amf, vol)
! dlnq(iclust) = emf / (kb*temp**2) +
! end associate
! end do
do iclust = 1, size(clusterset)
associate(c => clusterset(iclust))
emf = -real(sum(c%composition), dp)*sum(global_data%ntot)/vol * &
(amf - temp * amf_temp)
dlnq(iclust) = emf / (kb*temp**2)
end associate
end do
end subroutine calculate_dlnqint
!=================================================================================
! Calculates the second temperature derivative of the translational partition
! function.
subroutine calculate_ddlnqtrans(dlnq, bxv, bxv_temp, temp, vol)
real(dp), dimension(:), intent(out) :: dlnq
real(dp), intent(in) :: bxv
real(dp), intent(in) :: bxv_temp
real(dp), intent(in) :: temp
real(dp), intent(in) :: vol
dlnq(:) = -1.5_dp/temp**2 - &
(global_data%vexcl * bxv_temp)**2/(-vol+bxv*global_data%vexcl)**2
! write(*,*) temp, vol, bxv*global_data%vexcl
end subroutine calculate_ddlnqtrans
!=================================================================================
! Calculates the second temperature derivative of the vibrational partition
! function.
subroutine calculate_ddlnqvib(dlnq, temp)
real(dp), dimension(:), intent(out) :: dlnq
real(dp), intent(in) :: temp
integer :: iclust
integer:: ifreq
real(dp):: factor
real(dp):: t_vib
real(dp):: fsinh
real(dp):: fcosh
real(dp):: q_ho, q_fr
real(dp):: w
! TODO: Check atoms.
factor = planck*100.0_dp*speed_of_light/kb
do iclust = 1, size(clusterset)
associate(f => clusterset(iclust)%frequencies, q => dlnq(iclust), &
x => clusterset(iclust)%anharmonicity)
q = 0.0_dp
q_ho = 0.0_dp
q_fr = 0.0_dp
do ifreq = 1, size(f)
t_vib = factor*f(ifreq)
! Check if anharmonicity constant is set
if (x > 0.0_dp) then
! Morse oscillator.
fsinh = sinh(0.5_dp*t_vib/temp)
fcosh = cosh(0.5_dp*t_vib/temp)
q = q - fcosh/fsinh*t_vib/(temp**3) - &
2.0_dp*x*t_vib**2/(temp**4)*fcosh/fsinh**3 + &
(t_vib/(2.0_dp*temp**2*fsinh))**2 + &
2.0_dp*x*t_vib/(temp**3)*(0.25_dp + 0.5_dp/fsinh**2) + &
x*(t_vib**3)/(4.0_dp*temp**5)*(2.0_dp*fcosh**2+1.0_dp) / &
fsinh**4
else
! Harmonic oscillator.
! TODO: Use log1p and exp1m family of functions.
q_ho = -t_vib/(temp**3) &
- 2.0_dp*t_vib/(temp**3*(exp(t_vib/temp)-1.0_dp)) &
+ t_vib**2*exp(t_vib/temp) / &
(temp**4*(exp(t_vib/temp)-1)**2)
! Free rotator approximation by Grimme.
q_fr = -0.5_dp/temp**2
w = 1.0_dp/(1.0_dp + (global_data%rotor_cutoff/f(ifreq))**4)
q = q + w * q_ho + (1.0_dp - w) * q_fr
end if
end do
end associate
end do
end subroutine calculate_ddlnqvib
!=================================================================================
! Calculates the second temperature derivative of the rotational partition
! function.
subroutine calculate_ddlnqrot(dlnq, temp)
real(dp), dimension(:), intent(out) :: dlnq
real(dp), intent(in) :: temp
integer :: iclust
do iclust = 1, size(clusterset)
associate(c => clusterset(iclust), q => dlnq(iclust))
if (c%atom) then
q = 0.0_dp
else if (c%linear) then
! linear cluster:
! q_rot = T/(sigma*t_rot) -> d/dt(lnq)=-1/T**2
q = -1.0_dp/temp**2
else
! polyatomic, nonlinear cluster: d/dT(lnq)=-3/(2T**2)
q = -1.5_dp/temp**2
end if
end associate
end do
end subroutine calculate_ddlnqrot
!=================================================================================
! Calculates the second temperature derivative of the electronic partition
! function.
subroutine calculate_ddlnqelec(dlnq, temp)
real(dp), dimension(:), intent(out) :: dlnq
real(dp), intent(in) :: temp
! The adiabatic interaction energy is in units of kJ/mol.
dlnq(:) = - (2000.0_dp/avogadro)*clusterset(:)%energy/(kb*temp**3)
end subroutine calculate_ddlnqelec
!=================================================================================
! Calculates the second temperature derivative of the the mean field partition
! function.
subroutine calculate_ddlnqint(dlnq, temp, amf, amf_temp, vol)
real(dp), dimension(:), intent(out) :: dlnq
real(dp), intent(in) :: temp
real(dp), intent(in) :: amf
real(dp), intent(in) :: amf_temp
real(dp), intent(in) :: vol
integer :: iclust
real(dp):: emf
! do iclust = 1, size(clusterset)
! associate(c => clusterset(iclust))
! call calculate_emf(emf, c%composition, amf, vol)
! dlnq(iclust) = -2.0_dp*emf/(kb*temp**3)
! end associate
! end do
do iclust = 1, size(clusterset)
associate(c => clusterset(iclust))
emf = -real(sum(c%composition), dp)*sum(global_data%ntot)/vol * &
(amf - temp * amf_temp)
dlnq(iclust) = -2.0_dp*emf/(kb*temp**3)
end associate
end do
end subroutine calculate_ddlnqint
!=================================================================================
end module partition_functions
!=========================================================================================