Hide keyboard shortcuts

Hot-keys on this page

r m x p   toggle line displays

j k   next/prev highlighted chunk

0   (zero) top of page

1   (one) first highlighted chunk

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

1001

1002

1003

1004

1005

1006

1007

1008

1009

1010

1011

1012

1013

1014

1015

1016

1017

1018

1019

1020

1021

1022

1023

1024

1025

1026

1027

1028

1029

1030

1031

1032

1033

1034

1035

1036

1037

1038

1039

1040

1041

1042

1043

1044

1045

1046

1047

1048

1049

1050

1051

1052

1053

1054

1055

1056

1057

1058

1059

1060

1061

1062

1063

1064

1065

1066

1067

1068

1069

1070

1071

1072

1073

1074

1075

1076

1077

1078

1079

1080

1081

1082

1083

1084

1085

1086

1087

1088

1089

1090

1091

1092

1093

1094

1095

1096

1097

1098

1099

1100

1101

1102

1103

1104

1105

1106

1107

1108

1109

1110

1111

1112

1113

1114

1115

1116

1117

1118

1119

1120

1121

1122

1123

1124

1125

1126

1127

1128

1129

1130

1131

1132

1133

1134

1135

1136

1137

1138

1139

1140

1141

1142

1143

1144

1145

1146

1147

1148

1149

1150

1151

1152

1153

1154

1155

1156

1157

1158

1159

1160

1161

1162

1163

1164

1165

1166

1167

1168

1169

1170

1171

1172

1173

1174

1175

1176

1177

1178

1179

1180

1181

1182

1183

1184

1185

1186

1187

1188

1189

1190

1191

1192

1193

1194

1195

1196

1197

1198

1199

1200

1201

1202

1203

1204

1205

1206

1207

1208

1209

1210

1211

1212

1213

1214

1215

1216

1217

1218

1219

1220

1221

1222

1223

1224

1225

1226

1227

1228

1229

1230

1231

1232

1233

1234

1235

1236

1237

1238

1239

1240

1241

1242

1243

1244

1245

1246

1247

1248

1249

1250

1251

1252

1253

1254

1255

1256

1257

1258

1259

1260

1261

1262

1263

1264

1265

1266

1267

1268

1269

1270

1271

1272

1273

1274

1275

1276

1277

1278

1279

1280

1281

1282

1283

1284

1285

1286

1287

1288

1289

1290

1291

1292

1293

1294

1295

1296

1297

1298

1299

1300

1301

1302

1303

1304

1305

1306

1307

1308

1309

1310

1311

1312

1313

1314

1315

1316

1317

1318

1319

1320

1321

1322

1323

1324

1325

1326

1327

1328

1329

1330

1331

1332

1333

1334

1335

1336

1337

1338

1339

1340

1341

1342

1343

1344

1345

1346

1347

1348

1349

1350

1351

1352

1353

1354

1355

1356

1357

1358

1359

1360

1361

1362

1363

1364

1365

1366

1367

1368

1369

1370

1371

1372

1373

1374

1375

1376

1377

1378

1379

1380

1381

1382

1383

1384

1385

1386

1387

1388

1389

1390

1391

1392

1393

1394

1395

1396

1397

1398

1399

1400

1401

1402

1403

1404

1405

1406

1407

1408

1409

1410

1411

1412

1413

1414

1415

1416

1417

1418

1419

1420

1421

1422

1423

1424

1425

1426

1427

1428

1429

1430

1431

1432

1433

1434

1435

1436

1437

1438

1439

1440

1441

1442

1443

1444

1445

1446

1447

1448

1449

1450

1451

1452

1453

1454

1455

1456

# -*- coding: utf-8 -*- 

'''Chemical Engineering Design Library (ChEDL). Utilities for process modeling. 

Copyright (C) 2016, Caleb Bell <Caleb.Andrew.Bell@gmail.com> 

 

Permission is hereby granted, free of charge, to any person obtaining a copy 

of this software and associated documentation files (the "Software"), to deal 

in the Software without restriction, including without limitation the rights 

to use, copy, modify, merge, publish, distribute, sublicense, and/or sell 

copies of the Software, and to permit persons to whom the Software is 

furnished to do so, subject to the following conditions: 

 

The above copyright notice and this permission notice shall be included in all 

copies or substantial portions of the Software. 

 

THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR 

IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY, 

FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL THE 

AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER 

LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING FROM, 

OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN THE 

SOFTWARE.''' 

from __future__ import division 

 

__all__ = ['GCEOSMIX', 'PRMIX', 'SRKMIX', 'PR78MIX', 'VDWMIX', 'PRSVMIX', 

'PRSV2MIX', 'TWUPRMIX', 'TWUSRKMIX', 'APISRKMIX'] 

from scipy.constants import R 

from scipy.optimize import newton 

from thermo.utils import Cp_minus_Cv, isothermal_compressibility, phase_identification_parameter, phase_identification_parameter_phase 

from thermo.utils import log, exp, sqrt 

from thermo.utils import _isobaric_expansion as isobaric_expansion 

from thermo.eos import * 

 

 

class GCEOSMIX(GCEOS): 

r'''Class for solving a generic pressure-explicit three-parameter cubic  

equation of state for a mixture. Does not implement any parameters itself;  

must be subclassed by a mixture equation of state class which subclasses it. 

No routines for partial molar properties for a generic cubic equation of 

state have yet been implemented, although that would be desireable.  

The only partial molar property which is currently used is fugacity, which 

must be implemented in each mixture EOS that subclasses this. 

 

.. math:: 

P=\frac{RT}{V-b}-\frac{a\alpha(T)}{V^2 + \delta V + \epsilon} 

 

Main methods are `fugacities`, `solve_T`, and `a_alpha_and_derivatives`. 

 

`fugacities` is a helper method intended as a common interface for setting 

fugacities of each species in each phase; it calls `fugacity_coefficients` 

to actually calculate them, but that is not implemented here. This should 

be used when performing flash calculations, where fugacities are needed  

repeatedly. The fugacities change as a function of liquid/gas phase  

composition, but the entire EOS need not be solved to recalculate them. 

 

`solve_T` is a wrapper around `GCEOS`'s `solve_T`; the only difference is  

to use half the average mixture's critical temperature as the initial  

guess. 

 

`a_alpha_and_derivatives` implements the Van der Waals mixing rules for a 

mixture. It calls `a_alpha_and_derivatives` from the pure-component EOS for  

each species via multiple inheritance. 

''' 

def a_alpha_and_derivatives(self, T, full=True, quick=True): 

r'''Method to calculate `a_alpha` and its first and second 

derivatives for an EOS with the Van der Waals mixing rules. Uses the 

parent class's interface to compute pure component values. Returns 

`a_alpha`, `da_alpha_dT`, and `d2a_alpha_dT2`. Calls  

`setup_a_alpha_and_derivatives` before calling 

`a_alpha_and_derivatives` for each species, which typically sets `a`  

and `Tc`. Calls `cleanup_a_alpha_and_derivatives` to remove the set 

properties after the calls are done. 

 

For use in `solve_T` this returns only `a_alpha` if `full` is False. 

 

.. math:: 

a \alpha = \sum_i \sum_j z_i z_j {(a\alpha)}_{ij} 

 

(a\alpha)_{ij} = (1-k_{ij})\sqrt{(a\alpha)_{i}(a\alpha)_{j}} 

 

Parameters 

---------- 

T : float 

Temperature, [K] 

full : bool, optional 

If False, calculates and returns only `a_alpha` 

quick : bool, optional 

Only the quick variant is implemented; it is little faster anyhow 

 

Returns 

------- 

a_alpha : float 

Coefficient calculated by EOS-specific method, [J^2/mol^2/Pa] 

da_alpha_dT : float 

Temperature derivative of coefficient calculated by EOS-specific  

method, [J^2/mol^2/Pa/K] 

d2a_alpha_dT2 : float 

Second temperature derivative of coefficient calculated by  

EOS-specific method, [J^2/mol^2/Pa/K**2] 

 

Notes 

----- 

The exact expressions can be obtained with the following SymPy  

expression below, commented out for brevity. 

 

>>> from sympy import * 

>>> a_alpha_i, a_alpha_j, kij, T = symbols('a_alpha_i, a_alpha_j, kij, T') 

>>> a_alpha_ij = (1-kij)*sqrt(a_alpha_i(T)*a_alpha_j(T)) 

>>> #diff(a_alpha_ij, T) 

>>> #diff(a_alpha_ij, T, T) 

''' 

zs, kijs = self.zs, self.kijs 

a_alphas, da_alpha_dTs, d2a_alpha_dT2s = [], [], [] 

 

for i in self.cmps: 

self.setup_a_alpha_and_derivatives(i, T=T) 

# Abuse method resolution order to call the a_alpha_and_derivatives 

# method of the original pure EOS 

# -4 goes back from object, GCEOS, SINGLEPHASEEOS, up to GCEOSMIX 

ds = super(type(self).__mro__[self.a_alpha_mro], self).a_alpha_and_derivatives(T) 

a_alphas.append(ds[0]) 

da_alpha_dTs.append(ds[1]) 

d2a_alpha_dT2s.append(ds[2]) 

self.cleanup_a_alpha_and_derivatives() 

 

a_alpha, da_alpha_dT, d2a_alpha_dT2 = 0.0, 0.0, 0.0 

self.a_alpha_ijs = [[0]*self.N for i in self.cmps] 

 

for i in self.cmps: 

for j in self.cmps: 

self.a_alpha_ijs[i][j] = (1. - kijs[i][j])*(a_alphas[i]*a_alphas[j])**0.5 

# Needed in calculation of fugacity coefficients 

 

for i in self.cmps: 

for j in self.cmps: 

a_alpha += self.a_alpha_ijs[i][j]*zs[i]*zs[j] 

 

if full: 

for i in self.cmps: 

for j in self.cmps: 

da_alpha_dT += zs[i]*zs[j]*((1. - kijs[i][j])/(2*(a_alphas[i]*a_alphas[j])**0.5) 

*(a_alphas[i]*da_alpha_dTs[j] + a_alphas[j]*da_alpha_dTs[i])) 

 

x0 = a_alphas[i]*a_alphas[j] 

x1 = a_alphas[i]*da_alpha_dTs[j] 

x2 = a_alphas[j]*da_alpha_dTs[i] 

x3 = 2.*a_alphas[i]*da_alpha_dTs[j] + 2.*a_alphas[j]*da_alpha_dTs[i] 

d2a_alpha_dT2 += (-(x0)**0.5*(kijs[i][j] - 1.)*(x0*( 

2.*a_alphas[i]*d2a_alpha_dT2s[j] + 2.*a_alphas[j]*d2a_alpha_dT2s[i] 

+ 4.*da_alpha_dTs[i]*da_alpha_dTs[j]) - x1*x3 - x2*x3 + (x1 

+ x2)**2)/(4.*a_alphas[i]**2*a_alphas[j]**2))*zs[i]*zs[j] 

 

if full: 

return a_alpha, da_alpha_dT, d2a_alpha_dT2 

else: 

return a_alpha 

 

def fugacities(self, xs=None, ys=None): 

r'''Helper method for calculating fugacity coefficients for any  

phases present, using either the overall mole fractions for both phases 

or using specified mole fractions for each phase. 

 

Requires `fugacity_coefficients` to be implemented by each subclassing 

EOS. 

 

In addition to setting `fugacities_l` and/or `fugacities_g`, this also 

sets the fugacity coefficients `phis_l` and/or `phis_g`. 

 

.. math:: 

\hat \phi_i^g = \frac{\hat f_i^g}{x_i P} 

 

\hat \phi_i^l = \frac{\hat f_i^l}{x_i P} 

 

Parameters 

---------- 

xs : list[float], optional 

Liquid-phase mole fractions of each species, [-] 

ys : list[float], optional 

Vapor-phase mole fractions of each species, [-] 

 

Notes 

----- 

It is helpful to check that `fugacity_coefficients` has been 

implemented correctly using the following expression, from [1]_. 

 

.. math:: 

\ln \hat \phi_i = \left[\frac{\partial (n\log \phi)}{\partial  

n_i}\right]_{T,P,n_j,V_t} 

 

For reference, several expressions for fugacity of a component are as 

follows, shown in [1]_ and [2]_. 

 

.. math:: 

\ln \hat \phi_i = \int_{0}^P\left(\frac{\hat V_i} 

{RT} - \frac{1}{P}\right)dP 

 

\ln \hat \phi_i = \int_V^\infty \left[ 

\frac{1}{RT}\frac{\partial P}{ \partial n_i} 

- \frac{1}{V}\right] d V - \ln Z 

 

References 

---------- 

.. [1] Hu, Jiawen, Rong Wang, and Shide Mao. "Some Useful Expressions  

for Deriving Component Fugacity Coefficients from Mixture Fugacity  

Coefficient." Fluid Phase Equilibria 268, no. 1-2 (June 25, 2008):  

7-13. doi:10.1016/j.fluid.2008.03.007. 

.. [2] Walas, Stanley M. Phase Equilibria in Chemical Engineering.  

Butterworth-Heinemann, 1985. 

''' 

if self.phase in ['l', 'l/g']: 

if xs is None: 

xs = self.zs 

Z = self.P*self.V_l/(R*self.T) 

self.phis_l = self.fugacity_coefficients(Z, zs=xs) 

self.fugacities_l = [phi*x*self.P for phi, x in zip(self.phis_l, xs)] 

self.lnphis_l = [log(i) for i in self.phis_l] 

if self.phase in ['g', 'l/g']: 

if ys is None: 

ys = self.zs 

fs, phis = [], [] 

Z = self.P*self.V_g/(R*self.T) 

self.phis_g = self.fugacity_coefficients(Z, zs=ys) 

self.fugacities_g = [phi*y*self.P for phi, y in zip(self.phis_g, ys)] 

self.lnphis_g = [log(i) for i in self.phis_g] 

 

def solve_T(self, P, V, quick=True): 

r'''Generic method to calculate `T` from a specified `P` and `V`. 

Provides SciPy's `newton` solver, and iterates to solve the general 

equation for `P`, recalculating `a_alpha` as a function of temperature 

using `a_alpha_and_derivatives` each iteration. 

 

Parameters 

---------- 

P : float 

Pressure, [Pa] 

V : float 

Molar volume, [m^3/mol] 

quick : bool, optional 

Unimplemented, although it may be possible to derive explicit  

expressions as done for many pure-component EOS 

 

Returns 

------- 

T : float 

Temperature, [K] 

''' 

self.Tc = sum(self.Tcs)/self.N 

# -4 goes back from object, GCEOS 

return super(type(self).__mro__[-3], self).solve_T(P=P, V=V, quick=quick) 

 

 

 

class PRMIX(GCEOSMIX, PR): 

r'''Class for solving the Peng-Robinson cubic equation of state for a  

mixture of any number of compounds. Subclasses `PR`. Solves the EOS on 

initialization and calculates fugacities for all components in all phases. 

 

The implemented method here is `fugacity_coefficients`, which implements 

the formula for fugacity coefficients in a mixture as given in [1]_. 

Two of `T`, `P`, and `V` are needed to solve the EOS. 

 

.. math:: 

P = \frac{RT}{v-b}-\frac{a\alpha(T)}{v(v+b)+b(v-b)} 

 

a \alpha = \sum_i \sum_j z_i z_j {(a\alpha)}_{ij} 

 

(a\alpha)_{ij} = (1-k_{ij})\sqrt{(a\alpha)_{i}(a\alpha)_{j}} 

 

b = \sum_i z_i b_i 

 

a_i=0.45724\frac{R^2T_{c,i}^2}{P_{c,i}} 

 

b_i=0.07780\frac{RT_{c,i}}{P_{c,i}} 

 

\alpha(T)_i=[1+\kappa_i(1-\sqrt{T_{r,i}})]^2 

 

\kappa_i=0.37464+1.54226\omega_i-0.26992\omega^2_i 

 

Parameters 

---------- 

Tcs : float 

Critical temperatures of all compounds, [K] 

Pcs : float 

Critical pressures of all compounds, [Pa] 

omegas : float 

Acentric factors of all compounds, [-] 

zs : float 

Overall mole fractions of all species, [-] 

kijs : list[list[float]], optional 

n*n size list of lists with binary interaction parameters for the 

Van der Waals mixing rules, default all 0 [-] 

T : float, optional 

Temperature, [K] 

P : float, optional 

Pressure, [Pa] 

V : float, optional 

Molar volume, [m^3/mol] 

 

Examples 

-------- 

T-P initialization, nitrogen-methane at 115 K and 1 MPa: 

 

>>> eos = PRMIX(T=115, P=1E6, Tcs=[126.1, 190.6], Pcs=[33.94E5, 46.04E5], omegas=[0.04, 0.011], zs=[0.5, 0.5], kijs=[[0,0],[0,0]]) 

>>> eos.V_l, eos.V_g 

(3.625735065042031e-05, 0.0007006656856469095) 

>>> eos.fugacities_l, eos.fugacities_g 

([793860.8382114634, 73468.55225303846], [436530.9247009119, 358114.63827532396]) 

 

Notes 

----- 

For P-V initializations, SciPy's `newton` solver is used to find T. 

 

References 

---------- 

.. [1] Peng, Ding-Yu, and Donald B. Robinson. "A New Two-Constant Equation  

of State." Industrial & Engineering Chemistry Fundamentals 15, no. 1  

(February 1, 1976): 59-64. doi:10.1021/i160057a011. 

.. [2] Robinson, Donald B., Ding-Yu Peng, and Samuel Y-K Chung. "The  

Development of the Peng - Robinson Equation and Its Application to Phase 

Equilibrium in a System Containing Methanol." Fluid Phase Equilibria 24, 

no. 1 (January 1, 1985): 25-41. doi:10.1016/0378-3812(85)87035-7.  

''' 

a_alpha_mro = -4 

def __init__(self, Tcs, Pcs, omegas, zs, kijs=None, T=None, P=None, V=None): 

self.N = len(Tcs) 

self.cmps = range(self.N) 

self.Tcs = Tcs 

self.Pcs = Pcs 

self.omegas = omegas 

self.zs = zs 

if kijs is None: 

kijs = [[0]*self.N for i in range(self.N)] 

self.kijs = kijs 

self.T = T 

self.P = P 

self.V = V 

 

self.ais = [self.c1*R*R*Tc*Tc/Pc for Tc, Pc in zip(Tcs, Pcs)] 

self.bs = [self.c2*R*Tc/Pc for Tc, Pc in zip(Tcs, Pcs)] 

self.b = sum(bi*zi for bi, zi in zip(self.bs, self.zs)) 

self.kappas = [0.37464 + 1.54226*omega - 0.26992*omega*omega for omega in omegas] 

 

self.delta = 2.*self.b 

self.epsilon = -self.b*self.b 

 

self.solve() 

self.fugacities() 

 

def setup_a_alpha_and_derivatives(self, i, T=None): 

r'''Sets `a`, `kappa`, and `Tc` for a specific component before the  

pure-species EOS's `a_alpha_and_derivatives` method is called. Both are  

called by `GCEOSMIX.a_alpha_and_derivatives` for every component.''' 

self.a, self.kappa, self.Tc = self.ais[i], self.kappas[i], self.Tcs[i] 

 

def cleanup_a_alpha_and_derivatives(self): 

r'''Removes properties set by `setup_a_alpha_and_derivatives`; run by 

`GCEOSMIX.a_alpha_and_derivatives` after `a_alpha` is calculated for  

every component''' 

del(self.a, self.kappa, self.Tc) 

 

def fugacity_coefficients(self, Z, zs): 

r'''Literature formula for calculating fugacity coefficients for each 

species in a mixture. Verified numerically. Applicable to most  

derivatives of the Peng-Robinson equation of state as well. 

Called by `fugacities` on initialization, or by a solver routine  

which is performing a flash calculation. 

 

.. math:: 

\ln \hat \phi_i = \frac{B_i}{B}(Z-1)-\ln(Z-B) + \frac{A}{2\sqrt{2}B} 

\left[\frac{B_i}{B} - \frac{2}{a\alpha}\sum_i y_i(a\alpha)_{ij}\right] 

\log\left[\frac{Z + (1+\sqrt{2})B}{Z-(\sqrt{2}-1)B}\right] 

 

A = \frac{(a\alpha)P}{R^2 T^2} 

 

B = \frac{b P}{RT} 

 

Parameters 

---------- 

Z : float 

Compressibility of the mixture for a desired phase, [-] 

zs : list[float], optional 

List of mole factions, either overall or in a specific phase, [-] 

 

Returns 

------- 

phis : float 

Fugacity coefficient for each species, [-] 

 

References 

---------- 

.. [1] Peng, Ding-Yu, and Donald B. Robinson. "A New Two-Constant  

Equation of State." Industrial & Engineering Chemistry Fundamentals  

15, no. 1 (February 1, 1976): 59-64. doi:10.1021/i160057a011. 

.. [2] Walas, Stanley M. Phase Equilibria in Chemical Engineering.  

Butterworth-Heinemann, 1985. 

''' 

A = self.a_alpha*self.P/R**2/self.T**2 

B = self.b*self.P/R/self.T 

phis = [] 

for i in self.cmps: 

t1 = self.bs[i]/self.b*(Z-1) - log(Z-B) 

t2 = 2./self.a_alpha*sum([zs[j]*self.a_alpha_ijs[i][j] for j in self.cmps]) 

t3 = t1 - A/(2*2**0.5*B)*(t2 - self.bs[i]/self.b)*log((Z + (sqrt(2)+1)*B)/(Z - (sqrt(2)-1)*B)) 

phis.append(exp(t3)) 

return phis 

 

 

class SRKMIX(GCEOSMIX, SRK): 

r'''Class for solving the Soave-Redlich-Kwong cubic equation of state for a  

mixture of any number of compounds. Subclasses `SRK`. Solves the EOS on 

initialization and calculates fugacities for all components in all phases. 

 

The implemented method here is `fugacity_coefficients`, which implements 

the formula for fugacity coefficients in a mixture as given in [1]_. 

Two of `T`, `P`, and `V` are needed to solve the EOS. 

 

.. math:: 

P = \frac{RT}{V-b} - \frac{a\alpha(T)}{V(V+b)} 

 

a \alpha = \sum_i \sum_j z_i z_j {(a\alpha)}_{ij} 

 

(a\alpha)_{ij} = (1-k_{ij})\sqrt{(a\alpha)_{i}(a\alpha)_{j}} 

 

b = \sum_i z_i b_i 

 

a_i =\left(\frac{R^2(T_{c,i})^{2}}{9(\sqrt[3]{2}-1)P_{c,i}} \right) 

=\frac{0.42748\cdot R^2(T_{c,i})^{2}}{P_{c,i}} 

 

b_i =\left( \frac{(\sqrt[3]{2}-1)}{3}\right)\frac{RT_{c,i}}{P_{c,i}} 

=\frac{0.08664\cdot R T_{c,i}}{P_{c,i}} 

 

\alpha(T)_i = \left[1 + m_i\left(1 - \sqrt{\frac{T}{T_{c,i}}}\right)\right]^2 

 

m_i = 0.480 + 1.574\omega_i - 0.176\omega_i^2 

 

Parameters 

---------- 

Tcs : float 

Critical temperatures of all compounds, [K] 

Pcs : float 

Critical pressures of all compounds, [Pa] 

omegas : float 

Acentric factors of all compounds, [-] 

zs : float 

Overall mole fractions of all species, [-] 

kijs : list[list[float]], optional 

n*n size list of lists with binary interaction parameters for the 

Van der Waals mixing rules, default all 0 [-] 

T : float, optional 

Temperature, [K] 

P : float, optional 

Pressure, [Pa] 

V : float, optional 

Molar volume, [m^3/mol] 

 

Examples 

-------- 

T-P initialization, nitrogen-methane at 115 K and 1 MPa: 

 

>>> eos = SRKMIX(T=115, P=1E6, Tcs=[126.1, 190.6], Pcs=[33.94E5, 46.04E5], omegas=[0.04, 0.011], zs=[0.5, 0.5], kijs=[[0,0],[0,0]]) 

>>> eos.V_l, eos.V_g 

(4.104755570185178e-05, 0.0007110155639819184) 

>>> eos.fugacities_l, eos.fugacities_g 

([817841.6430546846, 72382.81925202628], [442137.1280124604, 361820.79211909405]) 

 

Notes 

----- 

For P-V initializations, SciPy's `newton` solver is used to find T. 

 

References 

---------- 

.. [1] Soave, Giorgio. "Equilibrium Constants from a Modified Redlich-Kwong 

Equation of State." Chemical Engineering Science 27, no. 6 (June 1972):  

1197-1203. doi:10.1016/0009-2509(72)80096-4. 

.. [2] Poling, Bruce E. The Properties of Gases and Liquids. 5th  

edition. New York: McGraw-Hill Professional, 2000. 

.. [3] Walas, Stanley M. Phase Equilibria in Chemical Engineering.  

Butterworth-Heinemann, 1985. 

''' 

a_alpha_mro = -4 

def __init__(self, Tcs, Pcs, omegas, zs, kijs=None, T=None, P=None, V=None): 

self.N = len(Tcs) 

self.cmps = range(self.N) 

self.Tcs = Tcs 

self.Pcs = Pcs 

self.omegas = omegas 

self.zs = zs 

if kijs is None: 

kijs = [[0]*self.N for i in range(self.N)] 

self.kijs = kijs 

self.T = T 

self.P = P 

self.V = V 

 

self.ais = [self.c1*R*R*Tc*Tc/Pc for Tc, Pc in zip(Tcs, Pcs)] 

self.bs = [self.c2*R*Tc/Pc for Tc, Pc in zip(Tcs, Pcs)] 

self.b = sum(bi*zi for bi, zi in zip(self.bs, self.zs)) 

self.ms = [0.480 + 1.574*omega - 0.176*omega*omega for omega in omegas] 

self.delta = self.b 

 

self.solve() 

self.fugacities() 

 

def setup_a_alpha_and_derivatives(self, i, T=None): 

r'''Sets `a`, `m`, and `Tc` for a specific component before the  

pure-species EOS's `a_alpha_and_derivatives` method is called. Both are  

called by `GCEOSMIX.a_alpha_and_derivatives` for every component.''' 

self.a, self.m, self.Tc = self.ais[i], self.ms[i], self.Tcs[i] 

 

def cleanup_a_alpha_and_derivatives(self): 

r'''Removes properties set by `setup_a_alpha_and_derivatives`; run by 

`GCEOSMIX.a_alpha_and_derivatives` after `a_alpha` is calculated for  

every component''' 

del(self.a, self.m, self.Tc) 

 

def fugacity_coefficients(self, Z, zs): 

r'''Literature formula for calculating fugacity coefficients for each 

species in a mixture. Verified numerically. Applicable to most  

derivatives of the SRK equation of state as well. 

Called by `fugacities` on initialization, or by a solver routine  

which is performing a flash calculation. 

 

.. math:: 

\ln \hat \phi_i = \frac{B_i}{B}(Z-1) - \ln(Z-B) + \frac{A}{B} 

\left[\frac{B_i}{B} - \frac{2}{a \alpha}\sum_i y_i(a\alpha)_{ij} 

\right]\ln\left(1+\frac{B}{Z}\right) 

 

A=\frac{a\alpha P}{R^2T^2} 

 

B = \frac{bP}{RT} 

 

Parameters 

---------- 

Z : float 

Compressibility of the mixture for a desired phase, [-] 

zs : list[float], optional 

List of mole factions, either overall or in a specific phase, [-] 

 

Returns 

------- 

phis : float 

Fugacity coefficient for each species, [-] 

 

References 

---------- 

.. [1] Soave, Giorgio. "Equilibrium Constants from a Modified  

Redlich-Kwong Equation of State." Chemical Engineering Science 27, 

no. 6 (June 1972): 1197-1203. doi:10.1016/0009-2509(72)80096-4. 

.. [2] Walas, Stanley M. Phase Equilibria in Chemical Engineering.  

Butterworth-Heinemann, 1985. 

''' 

A = self.a_alpha*self.P/R**2/self.T**2 

B = self.b*self.P/R/self.T 

phis = [] 

for i in self.cmps: 

Bi = self.bs[i]*self.P/R/self.T 

t1 = Bi/B*(Z-1) - log(Z - B) 

t2 = A/B*(Bi/B - 2./self.a_alpha*sum([zs[j]*self.a_alpha_ijs[i][j] for j in self.cmps])) 

t3 = log(1. + B/Z) 

t4 = t1 + t2*t3 

phis.append(exp(t4)) 

return phis 

 

 

class PR78MIX(PRMIX): 

r'''Class for solving the Peng-Robinson cubic equation of state for a  

mixture of any number of compounds according to the 1978 variant.  

Subclasses `PR`. Solves the EOS on initialization and calculates fugacities  

for all components in all phases. 

 

Two of `T`, `P`, and `V` are needed to solve the EOS. 

 

.. math:: 

P = \frac{RT}{v-b}-\frac{a\alpha(T)}{v(v+b)+b(v-b)} 

 

a \alpha = \sum_i \sum_j z_i z_j {(a\alpha)}_{ij} 

 

(a\alpha)_{ij} = (1-k_{ij})\sqrt{(a\alpha)_{i}(a\alpha)_{j}} 

 

b = \sum_i z_i b_i 

 

a_i=0.45724\frac{R^2T_{c,i}^2}{P_{c,i}} 

 

b_i=0.07780\frac{RT_{c,i}}{P_{c,i}} 

 

\alpha(T)_i=[1+\kappa_i(1-\sqrt{T_{r,i}})]^2 

 

\kappa_i = 0.37464+1.54226\omega_i-0.26992\omega_i^2 \text{ if } \omega_i 

\le 0.491 

 

\kappa_i = 0.379642 + 1.48503 \omega_i - 0.164423\omega_i^2 + 0.016666 

\omega_i^3 \text{ if } \omega_i > 0.491 

 

Parameters 

---------- 

Tcs : float 

Critical temperatures of all compounds, [K] 

Pcs : float 

Critical pressures of all compounds, [Pa] 

omegas : float 

Acentric factors of all compounds, [-] 

zs : float 

Overall mole fractions of all species, [-] 

kijs : list[list[float]], optional 

n*n size list of lists with binary interaction parameters for the 

Van der Waals mixing rules, default all 0 [-] 

T : float, optional 

Temperature, [K] 

P : float, optional 

Pressure, [Pa] 

V : float, optional 

Molar volume, [m^3/mol] 

 

Examples 

-------- 

T-P initialization, nitrogen-methane at 115 K and 1 MPa, with modified 

acentric factors to show the difference between `PRMIX` 

 

>>> eos = PR78MIX(T=115, P=1E6, Tcs=[126.1, 190.6], Pcs=[33.94E5, 46.04E5], omegas=[0.6, 0.7], zs=[0.5, 0.5], kijs=[[0,0],[0,0]]) 

>>> eos.V_l, eos.V_g 

(3.239642793468722e-05, 0.000504337849300222) 

>>> eos.fugacities_l, eos.fugacities_g 

([833048.4511980319, 6160.908815331634], [460717.2776793947, 279598.90103207633]) 

 

Notes 

----- 

This variant is recommended over the original. 

 

References 

---------- 

.. [1] Peng, Ding-Yu, and Donald B. Robinson. "A New Two-Constant Equation  

of State." Industrial & Engineering Chemistry Fundamentals 15, no. 1  

(February 1, 1976): 59-64. doi:10.1021/i160057a011. 

.. [2] Robinson, Donald B., Ding-Yu Peng, and Samuel Y-K Chung. "The  

Development of the Peng - Robinson Equation and Its Application to Phase 

Equilibrium in a System Containing Methanol." Fluid Phase Equilibria 24, 

no. 1 (January 1, 1985): 25-41. doi:10.1016/0378-3812(85)87035-7.  

''' 

a_alpha_mro = -4 

def __init__(self, Tcs, Pcs, omegas, zs, kijs=None, T=None, P=None, V=None): 

self.N = len(Tcs) 

self.cmps = range(self.N) 

self.Tcs = Tcs 

self.Pcs = Pcs 

self.omegas = omegas 

self.zs = zs 

if kijs is None: 

kijs = [[0]*self.N for i in range(self.N)] 

self.kijs = kijs 

self.T = T 

self.P = P 

self.V = V 

 

self.ais = [self.c1*R*R*Tc*Tc/Pc for Tc, Pc in zip(Tcs, Pcs)] 

self.bs = [self.c2*R*Tc/Pc for Tc, Pc in zip(Tcs, Pcs)] 

self.b = sum(bi*zi for bi, zi in zip(self.bs, self.zs)) 

self.kappas = [] 

for omega in omegas: 

if omega <= 0.491: 

self.kappas.append(0.37464 + 1.54226*omega - 0.26992*omega*omega) 

else: 

self.kappas.append(0.379642 + 1.48503*omega - 0.164423*omega**2 + 0.016666*omega**3) 

 

self.delta = 2.*self.b 

self.epsilon = -self.b*self.b 

 

self.solve() 

self.fugacities() 

 

 

 

class VDWMIX(GCEOSMIX, VDW): 

r'''Class for solving the Van der Waals cubic equation of state for a  

mixture of any number of compounds. Subclasses `VDW`. Solves the EOS on 

initialization and calculates fugacities for all components in all phases. 

 

The implemented method here is `fugacity_coefficients`, which implements 

the formula for fugacity coefficients in a mixture as given in [1]_. 

Two of `T`, `P`, and `V` are needed to solve the EOS. 

 

.. math:: 

P=\frac{RT}{V-b}-\frac{a}{V^2} 

 

a = \sum_i \sum_j z_i z_j {a}_{ij} 

 

b = \sum_i z_i b_i 

 

a_{ij} = (1-k_{ij})\sqrt{a_{i}a_{j}} 

 

a_i=\frac{27}{64}\frac{(RT_{c,i})^2}{P_{c,i}} 

 

b_i=\frac{RT_{c,i}}{8P_{c,i}} 

 

Parameters 

---------- 

Tcs : float 

Critical temperatures of all compounds, [K] 

Pcs : float 

Critical pressures of all compounds, [Pa] 

zs : float 

Overall mole fractions of all species, [-] 

kijs : list[list[float]], optional 

n*n size list of lists with binary interaction parameters for the 

Van der Waals mixing rules, default all 0 [-] 

T : float, optional 

Temperature, [K] 

P : float, optional 

Pressure, [Pa] 

V : float, optional 

Molar volume, [m^3/mol] 

 

Examples 

-------- 

T-P initialization, nitrogen-methane at 115 K and 1 MPa: 

 

>>> eos = VDWMIX(T=115, P=1E6, Tcs=[126.1, 190.6], Pcs=[33.94E5, 46.04E5], zs=[0.5, 0.5], kijs=[[0,0],[0,0]]) 

>>> eos.V_l, eos.V_g 

(5.8813678514166464e-05, 0.0007770869741895237) 

>>> eos.fugacities_l, eos.fugacities_g 

([854533.2669205102, 207126.8497276205], [448470.7363380735, 397826.5439999289]) 

 

Notes 

----- 

For P-V initializations, SciPy's `newton` solver is used to find T. 

 

References 

---------- 

.. [1] Walas, Stanley M. Phase Equilibria in Chemical Engineering.  

Butterworth-Heinemann, 1985. 

.. [2] Poling, Bruce E. The Properties of Gases and Liquids. 5th  

edition. New York: McGraw-Hill Professional, 2000. 

''' 

a_alpha_mro = -4 

def __init__(self, Tcs, Pcs, zs, kijs=None, T=None, P=None, V=None): 

self.N = len(Tcs) 

self.cmps = range(self.N) 

self.Tcs = Tcs 

self.Pcs = Pcs 

self.zs = zs 

if kijs is None: 

kijs = [[0]*self.N for i in range(self.N)] 

self.kijs = kijs 

self.T = T 

self.P = P 

self.V = V 

 

self.ais = [27.0/64.0*(R*Tc)**2/Pc for Tc, Pc in zip(Tcs, Pcs)] 

self.bs = [R*Tc/(8.*Pc) for Tc, Pc in zip(Tcs, Pcs)] 

self.b = sum(bi*zi for bi, zi in zip(self.bs, self.zs)) 

 

self.solve() 

self.fugacities() 

 

def setup_a_alpha_and_derivatives(self, i, T=None): 

r'''Sets `a` for a specific component before the  

pure-species EOS's `a_alpha_and_derivatives` method is called. Both are  

called by `GCEOSMIX.a_alpha_and_derivatives` for every component.''' 

self.a = self.ais[i] 

 

def cleanup_a_alpha_and_derivatives(self): 

r'''Removes properties set by `setup_a_alpha_and_derivatives`; run by 

`GCEOSMIX.a_alpha_and_derivatives` after `a_alpha` is calculated for  

every component''' 

del(self.a) 

 

def fugacity_coefficients(self, Z, zs): 

r'''Literature formula for calculating fugacity coefficients for each 

species in a mixture. Verified numerically. 

Called by `fugacities` on initialization, or by a solver routine  

which is performing a flash calculation. 

 

.. math:: 

\ln \hat \phi_i = \frac{b_i}{V-b} - \ln\left[Z\left(1 

- \frac{b}{V}\right)\right] - \frac{2\sqrt{aa_i}}{RTV} 

 

Parameters 

---------- 

Z : float 

Compressibility of the mixture for a desired phase, [-] 

zs : list[float], optional 

List of mole factions, either overall or in a specific phase, [-] 

 

Returns 

------- 

phis : float 

Fugacity coefficient for each species, [-] 

 

References 

---------- 

.. [1] Walas, Stanley M. Phase Equilibria in Chemical Engineering.  

Butterworth-Heinemann, 1985. 

''' 

phis = [] 

V = Z*R*self.T/self.P 

for i in self.cmps: 

phi = self.bs[i]/(V-self.b) - log(Z*(1. - self.b/V)) - 2.*(self.a_alpha*self.ais[i])**0.5/(R*self.T*V) 

phis.append(exp(phi)) 

return phis 

 

 

class PRSVMIX(PRMIX, PRSV): 

r'''Class for solving the Peng-Robinson-Stryjek-Vera equations of state for 

a mixture as given in [1]_. Subclasses `PRMIX` and `PRSV`. 

Solves the EOS on initialization and calculates fugacities for all  

components in all phases. 

 

Inherits the method of calculating fugacity coefficients from `PRMIX`. 

Two of `T`, `P`, and `V` are needed to solve the EOS. 

 

.. math:: 

P = \frac{RT}{v-b}-\frac{a\alpha(T)}{v(v+b)+b(v-b)} 

 

a \alpha = \sum_i \sum_j z_i z_j {(a\alpha)}_{ij} 

 

(a\alpha)_{ij} = (1-k_{ij})\sqrt{(a\alpha)_{i}(a\alpha)_{j}} 

 

b = \sum_i z_i b_i 

 

a_i=0.45724\frac{R^2T_{c,i}^2}{P_{c,i}} 

 

b_i=0.07780\frac{RT_{c,i}}{P_{c,i}} 

 

\alpha(T)_i=[1+\kappa_i(1-\sqrt{T_{r,i}})]^2 

 

\kappa_i = \kappa_{0,i} + \kappa_{1,i}(1 + T_{r,i}^{0.5})(0.7 - T_{r,i}) 

 

\kappa_{0,i} = 0.378893 + 1.4897153\omega_i - 0.17131848\omega_i^2  

+ 0.0196554\omega_i^3 

 

Parameters 

---------- 

Tcs : float 

Critical temperatures of all compounds, [K] 

Pcs : float 

Critical pressures of all compounds, [Pa] 

omegas : float 

Acentric factors of all compounds, [-] 

zs : float 

Overall mole fractions of all species, [-] 

kijs : list[list[float]], optional 

n*n size list of lists with binary interaction parameters for the 

Van der Waals mixing rules, default all 0 [-] 

T : float, optional 

Temperature, [K] 

P : float, optional 

Pressure, [Pa] 

V : float, optional 

Molar volume, [m^3/mol] 

kappa1s : list[float], optional 

Fit parameter; available in [1]_ for over 90 compounds, [-] 

 

Examples 

-------- 

P-T initialization, two-phase, nitrogen and methane 

 

>>> eos = PRSVMIX(T=115, P=1E6, Tcs=[126.1, 190.6], Pcs=[33.94E5, 46.04E5], omegas=[0.04, 0.011], zs=[0.5, 0.5], kijs=[[0,0],[0,0]]) 

>>> eos.phase, eos.V_l, eos.H_dep_l, eos.S_dep_l 

('l/g', 3.623552388375633e-05, -6349.003406339961, -49.12403359687138) 

 

Notes 

----- 

[1]_ recommends that `kappa1` be set to 0 for Tr > 0.7. This is not done by  

default; the class boolean `kappa1_Tr_limit` may be set to True and the 

problem re-solved with that specified if desired. `kappa1_Tr_limit` is not 

supported for P-V inputs. 

 

For P-V initializations, SciPy's `newton` solver is used to find T. 

 

[2]_ and [3]_ are two more resources documenting the PRSV EOS. [4]_ lists 

`kappa` values for 69 additional compounds. See also `PRSV2`. Note that 

tabulated `kappa` values should be used with the critical parameters used 

in their fits. Both [1]_ and [4]_ only considered vapor pressure in fitting 

the parameter. 

 

References 

---------- 

.. [1] Stryjek, R., and J. H. Vera. "PRSV: An Improved Peng-Robinson  

Equation of State for Pure Compounds and Mixtures." The Canadian Journal 

of Chemical Engineering 64, no. 2 (April 1, 1986): 323-33.  

doi:10.1002/cjce.5450640224.  

.. [2] Stryjek, R., and J. H. Vera. "PRSV - An Improved Peng-Robinson  

Equation of State with New Mixing Rules for Strongly Nonideal Mixtures." 

The Canadian Journal of Chemical Engineering 64, no. 2 (April 1, 1986):  

334-40. doi:10.1002/cjce.5450640225.  

.. [3] Stryjek, R., and J. H. Vera. "Vapor-liquid Equilibrium of  

Hydrochloric Acid Solutions with the PRSV Equation of State." Fluid  

Phase Equilibria 25, no. 3 (January 1, 1986): 279-90.  

doi:10.1016/0378-3812(86)80004-8.  

.. [4] Proust, P., and J. H. Vera. "PRSV: The Stryjek-Vera Modification of  

the Peng-Robinson Equation of State. Parameters for Other Pure Compounds 

of Industrial Interest." The Canadian Journal of Chemical Engineering  

67, no. 1 (February 1, 1989): 170-73. doi:10.1002/cjce.5450670125. 

''' 

a_alpha_mro = -5 

def __init__(self, Tcs, Pcs, omegas, zs, kijs=None, T=None, P=None, V=None, kappa1s=None): 

self.N = len(Tcs) 

self.cmps = range(self.N) 

self.Tcs = Tcs 

self.Pcs = Pcs 

self.omegas = omegas 

self.zs = zs 

 

if kijs is None: 

kijs = [[0]*self.N for i in range(self.N)] 

self.kijs = kijs 

 

if kappa1s is None: 

kappa1s = [0 for i in self.cmps] 

 

self.T = T 

self.P = P 

self.V = V 

 

self.ais = [self.c1*R*R*Tc*Tc/Pc for Tc, Pc in zip(Tcs, Pcs)] 

self.bs = [self.c2*R*Tc/Pc for Tc, Pc in zip(Tcs, Pcs)] 

self.b = sum(bi*zi for bi, zi in zip(self.bs, self.zs)) 

 

self.kappa0s = [0.378893 + 1.4897153*omega - 0.17131848*omega**2 + 0.0196554*omega**3 for omega in omegas] 

 

self.delta = 2*self.b 

self.epsilon = -self.b*self.b 

 

self.check_sufficient_inputs() 

if self.V and self.P: 

# Deal with T-solution here; does NOT support kappa1_Tr_limit. 

self.kappa1s = kappa1s 

self.T = self.solve_T(self.P, self.V) 

else: 

self.kappa1s = [(0 if (T/Tc > 0.7 and self.kappa1_Tr_limit) else kappa1) for kappa1, Tc in zip(kappa1s, Tcs)] 

 

self.kappas = [kappa0 + kappa1*(1 + (self.T/Tc)**0.5)*(0.7 - (self.T/Tc)) for kappa0, kappa1, Tc in zip(self.kappa0s, self.kappa1s, self.Tcs)] 

self.solve() 

 

self.fugacities() 

 

def setup_a_alpha_and_derivatives(self, i, T=None): 

r'''Sets `a`, `kappa0`, `kappa1`, and `Tc` for a specific component before the  

pure-species EOS's `a_alpha_and_derivatives` method is called. Both are  

called by `GCEOSMIX.a_alpha_and_derivatives` for every component.''' 

if not hasattr(self, 'kappas'): 

self.kappas = [kappa0 + kappa1*(1 + (T/Tc)**0.5)*(0.7 - (T/Tc)) for kappa0, kappa1, Tc in zip(self.kappa0s, self.kappa1s, self.Tcs)] 

self.a, self.kappa, self.kappa0, self.kappa1, self.Tc = self.ais[i], self.kappas[i], self.kappa0s[i], self.kappa1s[i], self.Tcs[i] 

 

def cleanup_a_alpha_and_derivatives(self): 

r'''Removes properties set by `setup_a_alpha_and_derivatives`; run by 

`GCEOSMIX.a_alpha_and_derivatives` after `a_alpha` is calculated for  

every component''' 

del(self.a, self.kappa, self.kappa0, self.kappa1, self.Tc) 

 

 

class PRSV2MIX(PRMIX, PRSV2): 

r'''Class for solving the Peng-Robinson-Stryjek-Vera 2 equations of state  

for a Mixture as given in [1]_. Subclasses `PRMIX` and `PRSV2`. 

Solves the EOS on initialization and calculates fugacities for all  

components in all phases. 

 

Inherits the method of calculating fugacity coefficients from `PRMIX`. 

Two of `T`, `P`, and `V` are needed to solve the EOS. 

 

.. math:: 

P = \frac{RT}{v-b}-\frac{a\alpha(T)}{v(v+b)+b(v-b)} 

 

a \alpha = \sum_i \sum_j z_i z_j {(a\alpha)}_{ij} 

 

(a\alpha)_{ij} = (1-k_{ij})\sqrt{(a\alpha)_{i}(a\alpha)_{j}} 

 

b = \sum_i z_i b_i 

 

a_i=0.45724\frac{R^2T_{c,i}^2}{P_{c,i}} 

 

b_i=0.07780\frac{RT_{c,i}}{P_{c,i}} 

 

\alpha(T)_i=[1+\kappa_i(1-\sqrt{T_{r,i}})]^2 

 

\kappa_i = \kappa_{0,i} + [\kappa_{1,i} + \kappa_{2,i}(\kappa_{3,i} - T_{r,i})(1-T_{r,i}^{0.5})] 

(1 + T_{r,i}^{0.5})(0.7 - T_{r,i}) 

 

\kappa_{0,i} = 0.378893 + 1.4897153\omega_i - 0.17131848\omega_i^2  

+ 0.0196554\omega_i^3 

 

Parameters 

---------- 

Tcs : float 

Critical temperatures of all compounds, [K] 

Pcs : float 

Critical pressures of all compounds, [Pa] 

omegas : float 

Acentric factors of all compounds, [-] 

zs : float 

Overall mole fractions of all species, [-] 

kijs : list[list[float]], optional 

n*n size list of lists with binary interaction parameters for the 

Van der Waals mixing rules, default all 0 [-] 

T : float, optional 

Temperature, [K] 

P : float, optional 

Pressure, [Pa] 

V : float, optional 

Molar volume, [m^3/mol] 

kappa1s : list[float], optional 

Fit parameter; available in [1]_ for over 90 compounds, [-] 

kappa2s : list[float], optional 

Fit parameter; available in [1]_ for over 90 compounds, [-] 

kappa3s : list[float], optional 

Fit parameter; available in [1]_ for over 90 compounds, [-] 

 

Examples 

-------- 

T-P initialization, nitrogen-methane at 115 K and 1 MPa: 

 

>>> eos = PRSV2MIX(T=115, P=1E6, Tcs=[126.1, 190.6], Pcs=[33.94E5, 46.04E5], omegas=[0.04, 0.011], zs=[0.5, 0.5], kijs=[[0,0],[0,0]]) 

>>> eos.V_l, eos.V_g 

(3.623552388375633e-05, 0.0007002421492037557) 

>>> eos.fugacities_l, eos.fugacities_g 

([794057.5831840546, 72851.22327178407], [436553.6561835047, 357878.1106688996]) 

 

Notes 

-----  

For P-V initializations, SciPy's `newton` solver is used to find T. 

 

Note that tabulated `kappa` values should be used with the critical  

parameters used in their fits. [1]_ considered only vapor  

pressure in fitting the parameter. 

 

References 

---------- 

.. [1] Stryjek, R., and J. H. Vera. "PRSV2: A Cubic Equation of State for  

Accurate Vapor-liquid Equilibria Calculations." The Canadian Journal of  

Chemical Engineering 64, no. 5 (October 1, 1986): 820-26.  

doi:10.1002/cjce.5450640516.  

''' 

a_alpha_mro = -5 

def __init__(self, Tcs, Pcs, omegas, zs, kijs=None, T=None, P=None, V=None, 

kappa1s=None, kappa2s=None, kappa3s=None): 

self.N = len(Tcs) 

self.cmps = range(self.N) 

self.Tcs = Tcs 

self.Pcs = Pcs 

self.omegas = omegas 

self.zs = zs 

 

if kijs is None: 

kijs = [[0]*self.N for i in range(self.N)] 

self.kijs = kijs 

 

if kappa1s is None: 

kappa1s = [0 for i in self.cmps] 

if kappa2s is None: 

kappa2s = [0 for i in self.cmps] 

if kappa3s is None: 

kappa3s = [0 for i in self.cmps] 

 

self.kappa1s = kappa1s 

self.kappa2s = kappa2s 

self.kappa3s = kappa3s 

 

self.T = T 

self.P = P 

self.V = V 

 

self.ais = [self.c1*R*R*Tc*Tc/Pc for Tc, Pc in zip(Tcs, Pcs)] 

self.bs = [self.c2*R*Tc/Pc for Tc, Pc in zip(Tcs, Pcs)] 

self.b = sum(bi*zi for bi, zi in zip(self.bs, self.zs)) 

 

self.kappa0s = [0.378893 + 1.4897153*omega - 0.17131848*omega**2 + 0.0196554*omega**3 for omega in omegas] 

 

self.delta = 2*self.b 

self.epsilon = -self.b*self.b 

 

 

if self.V and self.P: 

self.T = self.solve_T(self.P, self.V) 

 

self.kappas = [] 

for Tc, kappa0, kappa1, kappa2, kappa3 in zip(Tcs, self.kappa0s, self.kappa1s, self.kappa2s, self.kappa3s): 

Tr = self.T/Tc 

kappa = kappa0 + ((kappa1 + kappa2*(kappa3 - Tr)*(1. - Tr**0.5))*(1. + Tr**0.5)*(0.7 - Tr)) 

self.kappas.append(kappa) 

self.solve() 

self.fugacities() 

 

def setup_a_alpha_and_derivatives(self, i, T=None): 

r'''Sets `a`, `kappa`, `kappa0`, `kappa1`, `kappa2`, `kappa3` and `Tc` 

for a specific component before the  

pure-species EOS's `a_alpha_and_derivatives` method is called. Both are  

called by `GCEOSMIX.a_alpha_and_derivatives` for every component.''' 

if not hasattr(self, 'kappas'): 

self.kappas = [] 

for Tc, kappa0, kappa1, kappa2, kappa3 in zip(self.Tcs, self.kappa0s, self.kappa1s, self.kappa2s, self.kappa3s): 

Tr = T/Tc 

kappa = kappa0 + ((kappa1 + kappa2*(kappa3 - Tr)*(1. - Tr**0.5))*(1. + Tr**0.5)*(0.7 - Tr)) 

self.kappas.append(kappa) 

 

(self.a, self.kappa, self.kappa0, self.kappa1, self.kappa2, 

self.kappa3, self.Tc) = (self.ais[i], self.kappas[i], self.kappa0s[i], 

self.kappa1s[i], self.kappa2s[i], self.kappa3s[i], self.Tcs[i]) 

 

def cleanup_a_alpha_and_derivatives(self): 

r'''Removes properties set by `setup_a_alpha_and_derivatives`; run by 

`GCEOSMIX.a_alpha_and_derivatives` after `a_alpha` is calculated for  

every component''' 

del(self.a, self.kappa, self.kappa0, self.kappa1, self.kappa2, self.kappa3, self.Tc) 

 

 

class TWUPRMIX(PRMIX, TWUPR): 

r'''Class for solving the Twu [1]_ variant of the Peng-Robinson cubic  

equation of state for a mixture. Subclasses `TWUPR`. Solves the EOS on 

initialization and calculates fugacities for all components in all phases. 

 

Two of `T`, `P`, and `V` are needed to solve the EOS. 

 

.. math:: 

P = \frac{RT}{v-b}-\frac{a\alpha(T)}{v(v+b)+b(v-b)} 

 

a \alpha = \sum_i \sum_j z_i z_j {(a\alpha)}_{ij} 

 

(a\alpha)_{ij} = (1-k_{ij})\sqrt{(a\alpha)_{i}(a\alpha)_{j}} 

 

b = \sum_i z_i b_i 

 

a_i=0.45724\frac{R^2T_{c,i}^2}{P_{c,i}} 

 

b_i=0.07780\frac{RT_{c,i}}{P_{c,i}} 

 

\alpha_i = \alpha_i^{(0)} + \omega_i(\alpha_i^{(1)}-\alpha_i^{(0)}) 

 

\alpha^{(\text{0 or 1})} = T_{r,i}^{N(M-1)}\exp[L(1-T_{r,i}^{NM})] 

 

For sub-critical conditions: 

 

L0, M0, N0 = 0.125283, 0.911807, 1.948150; 

 

L1, M1, N1 = 0.511614, 0.784054, 2.812520 

 

For supercritical conditions: 

 

L0, M0, N0 = 0.401219, 4.963070, -0.2; 

 

L1, M1, N1 = 0.024955, 1.248089, -8.  

 

Parameters 

---------- 

Tcs : float 

Critical temperatures of all compounds, [K] 

Pcs : float 

Critical pressures of all compounds, [Pa] 

omegas : float 

Acentric factors of all compounds, [-] 

zs : float 

Overall mole fractions of all species, [-] 

kijs : list[list[float]], optional 

n*n size list of lists with binary interaction parameters for the 

Van der Waals mixing rules, default all 0 [-] 

T : float, optional 

Temperature, [K] 

P : float, optional 

Pressure, [Pa] 

V : float, optional 

Molar volume, [m^3/mol] 

 

Examples 

-------- 

T-P initialization, nitrogen-methane at 115 K and 1 MPa: 

 

>>> eos = TWUPRMIX(T=115, P=1E6, Tcs=[126.1, 190.6], Pcs=[33.94E5, 46.04E5], omegas=[0.04, 0.011], zs=[0.5, 0.5], kijs=[[0,0],[0,0]]) 

>>> eos.V_l, eos.V_g 

(3.62456981315702e-05, 0.0007004398944116554) 

>>> eos.fugacities_l, eos.fugacities_g 

([792155.0221633187, 73305.88829726784], [436468.96776424424, 358049.2495573095]) 

 

Notes 

----- 

For P-V initializations, SciPy's `newton` solver is used to find T. 

Claimed to be more accurate than the PR, PR78 and PRSV equations. 

 

References 

---------- 

.. [1] Twu, Chorng H., John E. Coon, and John R. Cunningham. "A New  

Generalized Alpha Function for a Cubic Equation of State Part 1.  

Peng-Robinson Equation." Fluid Phase Equilibria 105, no. 1 (March 15,  

1995): 49-59. doi:10.1016/0378-3812(94)02601-V. 

''' 

a_alpha_mro = -5 

def __init__(self, Tcs, Pcs, omegas, zs, kijs=None, T=None, P=None, V=None): 

self.N = len(Tcs) 

self.cmps = range(self.N) 

self.Tcs = Tcs 

self.Pcs = Pcs 

self.omegas = omegas 

self.zs = zs 

if kijs is None: 

kijs = [[0]*self.N for i in range(self.N)] 

self.kijs = kijs 

self.T = T 

self.P = P 

self.V = V 

 

self.ais = [self.c1*R*R*Tc*Tc/Pc for Tc, Pc in zip(Tcs, Pcs)] 

self.bs = [self.c2*R*Tc/Pc for Tc, Pc in zip(Tcs, Pcs)] 

self.b = sum(bi*zi for bi, zi in zip(self.bs, self.zs)) 

 

self.delta = 2.*self.b 

self.epsilon = -self.b*self.b 

self.check_sufficient_inputs() 

 

self.solve() 

self.fugacities() 

 

def setup_a_alpha_and_derivatives(self, i, T=None): 

r'''Sets `a`, `omega`, and `Tc` for a specific component before the  

pure-species EOS's `a_alpha_and_derivatives` method is called. Both are  

called by `GCEOSMIX.a_alpha_and_derivatives` for every component.''' 

self.a, self.Tc, self.omega = self.ais[i], self.Tcs[i], self.omegas[i] 

def cleanup_a_alpha_and_derivatives(self): 

r'''Removes properties set by `setup_a_alpha_and_derivatives`; run by 

`GCEOSMIX.a_alpha_and_derivatives` after `a_alpha` is calculated for  

every component''' 

del(self.a, self.Tc, self.omega) 

 

 

class TWUSRKMIX(SRKMIX, TWUSRK): 

r'''Class for solving the Twu variant of the Soave-Redlich-Kwong cubic  

equation of state for a mixture. Subclasses `TWUSRK`. Solves the EOS on 

initialization and calculates fugacities for all components in all phases. 

 

Two of `T`, `P`, and `V` are needed to solve the EOS. 

 

.. math:: 

P = \frac{RT}{V-b} - \frac{a\alpha(T)}{V(V+b)} 

 

a_i =\left(\frac{R^2(T_{c,i})^{2}}{9(\sqrt[3]{2}-1)P_{c,i}} \right) 

=\frac{0.42748\cdot R^2(T_{c,i})^{2}}{P_{c,i}} 

 

b_i =\left( \frac{(\sqrt[3]{2}-1)}{3}\right)\frac{RT_{c,i}}{P_{c,i}} 

=\frac{0.08664\cdot R T_{c,i}}{P_{c,i}} 

 

a \alpha = \sum_i \sum_j z_i z_j {(a\alpha)}_{ij} 

 

(a\alpha)_{ij} = (1-k_{ij})\sqrt{(a\alpha)_{i}(a\alpha)_{j}} 

 

b = \sum_i z_i b_i 

 

\alpha_i = \alpha^{(0,i)} + \omega_i(\alpha^{(1,i)}-\alpha^{(0,i)}) 

 

\alpha^{(\text{0 or 1, i})} = T_{r,i}^{N(M-1)}\exp[L(1-T_{r,i}^{NM})] 

 

For sub-critical conditions: 

 

L0, M0, N0 = 0.141599, 0.919422, 2.496441 

 

L1, M1, N1 = 0.500315, 0.799457, 3.291790 

 

For supercritical conditions: 

 

L0, M0, N0 = 0.441411, 6.500018, -0.20 

 

L1, M1, N1 = 0.032580, 1.289098, -8.0 

 

Parameters 

---------- 

Tcs : float 

Critical temperatures of all compounds, [K] 

Pcs : float 

Critical pressures of all compounds, [Pa] 

omegas : float 

Acentric factors of all compounds, [-] 

zs : float 

Overall mole fractions of all species, [-] 

kijs : list[list[float]], optional 

n*n size list of lists with binary interaction parameters for the 

Van der Waals mixing rules, default all 0 [-] 

T : float, optional 

Temperature, [K] 

P : float, optional 

Pressure, [Pa] 

V : float, optional 

Molar volume, [m^3/mol] 

 

Examples 

--------  

T-P initialization, nitrogen-methane at 115 K and 1 MPa: 

 

>>> eos = TWUSRKMIX(T=115, P=1E6, Tcs=[126.1, 190.6], Pcs=[33.94E5, 46.04E5], omegas=[0.04, 0.011], zs=[0.5, 0.5], kijs=[[0,0],[0,0]]) 

>>> eos.V_l, eos.V_g 

(4.108791361639091e-05, 0.0007117070840276789) 

>>> eos.fugacities_l, eos.fugacities_g 

([809692.8308266952, 74093.63881572781], [441783.43148985505, 362470.31741077645]) 

 

Notes 

----- 

For P-V initializations, SciPy's `newton` solver is used to find T. 

Claimed to be more accurate than the SRK equation. 

 

References 

---------- 

.. [1] Twu, Chorng H., John E. Coon, and John R. Cunningham. "A New  

Generalized Alpha Function for a Cubic Equation of State Part 2.  

Redlich-Kwong Equation." Fluid Phase Equilibria 105, no. 1 (March 15,  

1995): 61-69. doi:10.1016/0378-3812(94)02602-W. 

''' 

a_alpha_mro = -5 

def __init__(self, Tcs, Pcs, omegas, zs, kijs=None, T=None, P=None, V=None): 

self.N = len(Tcs) 

self.cmps = range(self.N) 

self.Tcs = Tcs 

self.Pcs = Pcs 

self.omegas = omegas 

self.zs = zs 

if kijs is None: 

kijs = [[0]*self.N for i in range(self.N)] 

self.kijs = kijs 

self.T = T 

self.P = P 

self.V = V 

 

self.ais = [self.c1*R*R*Tc*Tc/Pc for Tc, Pc in zip(Tcs, Pcs)] 

self.bs = [self.c2*R*Tc/Pc for Tc, Pc in zip(Tcs, Pcs)] 

self.b = sum(bi*zi for bi, zi in zip(self.bs, self.zs)) 

 

self.delta = self.b 

self.check_sufficient_inputs() 

 

self.solve() 

self.fugacities() 

 

def setup_a_alpha_and_derivatives(self, i, T=None): 

r'''Sets `a`, `omega`, and `Tc` for a specific component before the  

pure-species EOS's `a_alpha_and_derivatives` method is called. Both are  

called by `GCEOSMIX.a_alpha_and_derivatives` for every component.''' 

self.a, self.Tc, self.omega = self.ais[i], self.Tcs[i], self.omegas[i] 

 

def cleanup_a_alpha_and_derivatives(self): 

r'''Removes properties set by `setup_a_alpha_and_derivatives`; run by 

`GCEOSMIX.a_alpha_and_derivatives` after `a_alpha` is calculated for  

every component''' 

del(self.a, self.Tc, self.omega) 

 

 

class APISRKMIX(SRKMIX, APISRK): 

r'''Class for solving the Refinery Soave-Redlich-Kwong cubic  

equation of state for a mixture of any number of compounds, as shown in the 

API Databook [1]_. Subclasses `APISRK`. Solves the EOS on 

initialization and calculates fugacities for all components in all phases. 

 

Two of `T`, `P`, and `V` are needed to solve the EOS. 

 

.. math:: 

P = \frac{RT}{V-b} - \frac{a\alpha(T)}{V(V+b)} 

 

a \alpha = \sum_i \sum_j z_i z_j {(a\alpha)}_{ij} 

 

(a\alpha)_{ij} = (1-k_{ij})\sqrt{(a\alpha)_{i}(a\alpha)_{j}} 

 

b = \sum_i z_i b_i 

 

a_i =\left(\frac{R^2(T_{c,i})^{2}}{9(\sqrt[3]{2}-1)P_{c,i}} \right) 

=\frac{0.42748\cdot R^2(T_{c,i})^{2}}{P_{c,i}} 

 

b_i =\left( \frac{(\sqrt[3]{2}-1)}{3}\right)\frac{RT_{c,i}}{P_{c,i}} 

=\frac{0.08664\cdot R T_{c,i}}{P_{c,i}} 

 

\alpha(T)_i = \left[1 + S_{1,i}\left(1-\sqrt{T_{r,i}}\right) + S_{2,i} 

\frac{1- \sqrt{T_{r,i}}}{\sqrt{T_{r,i}}}\right]^2 

 

S_{1,i} = 0.48508 + 1.55171\omega_i - 0.15613\omega_i^2 \text{ if S1 is not tabulated } 

 

Parameters 

---------- 

Tcs : float 

Critical temperatures of all compounds, [K] 

Pcs : float 

Critical pressures of all compounds, [Pa] 

omegas : float 

Acentric factors of all compounds, [-] 

zs : float 

Overall mole fractions of all species, [-] 

kijs : list[list[float]], optional 

n*n size list of lists with binary interaction parameters for the 

Van der Waals mixing rules, default all 0 [-] 

T : float, optional 

Temperature, [K] 

P : float, optional 

Pressure, [Pa] 

V : float, optional 

Molar volume, [m^3/mol] 

S1s : float, optional 

Fit constant or estimated from acentric factor if not provided [-] 

S2s : float, optional 

Fit constant or 0 if not provided [-] 

 

Notes 

----- 

For P-V initializations, SciPy's `newton` solver is used to find T. 

 

Examples 

--------  

T-P initialization, nitrogen-methane at 115 K and 1 MPa: 

 

>>> eos = APISRKMIX(T=115, P=1E6, Tcs=[126.1, 190.6], Pcs=[33.94E5, 46.04E5], omegas=[0.04, 0.011], zs=[0.5, 0.5], kijs=[[0,0],[0,0]]) 

>>> eos.V_l, eos.V_g 

(4.101590920556748e-05, 0.0007104685894929316) 

>>> eos.fugacities_l, eos.fugacities_g 

([817882.3033490349, 71620.48238123364], [442158.29113191745, 361519.7987757053]) 

 

References 

---------- 

.. [1] API Technical Data Book: General Properties & Characterization. 

American Petroleum Institute, 7E, 2005. 

''' 

a_alpha_mro = -5 

def __init__(self, Tcs, Pcs, zs, omegas=None, kijs=None, T=None, P=None, V=None, 

S1s=None, S2s=None): 

self.N = len(Tcs) 

self.cmps = range(self.N) 

self.Tcs = Tcs 

self.Pcs = Pcs 

self.omegas = omegas 

self.zs = zs 

if kijs is None: 

kijs = [[0]*self.N for i in range(self.N)] 

self.kijs = kijs 

self.T = T 

self.P = P 

self.V = V 

 

self.check_sufficient_inputs() 

 

# Setup S1s and S2s 

if S1s is None and omegas is None: 

raise Exception('Either acentric factor of S1 is required') 

if S1s is None: 

self.S1s = [0.48508 + 1.55171*omega - 0.15613*omega*omega for omega in omegas] 

else: 

self.S1s = S1s 

if S2s is None: 

S2s = [0 for i in self.cmps] 

self.S2s = S2s 

 

self.ais = [self.c1*R*R*Tc*Tc/Pc for Tc, Pc in zip(Tcs, Pcs)] 

self.bs = [self.c2*R*Tc/Pc for Tc, Pc in zip(Tcs, Pcs)] 

self.b = sum(bi*zi for bi, zi in zip(self.bs, self.zs)) 

self.delta = self.b 

 

self.solve() 

self.fugacities() 

 

def setup_a_alpha_and_derivatives(self, i, T=None): 

r'''Sets `a`, `S1`, `S2` and `Tc` for a specific component before the  

pure-species EOS's `a_alpha_and_derivatives` method is called. Both are  

called by `GCEOSMIX.a_alpha_and_derivatives` for every component.''' 

self.a, self.Tc, self.S1, self.S2 = self.ais[i], self.Tcs[i], self.S1s[i], self.S2s[i] 

 

def cleanup_a_alpha_and_derivatives(self): 

r'''Removes properties set by `setup_a_alpha_and_derivatives`; run by 

`GCEOSMIX.a_alpha_and_derivatives` after `a_alpha` is calculated for  

every component''' 

del(self.a, self.Tc, self.S1, self.S2)