polykin.copolymerization¤
ImplicitPenultimateModel ¤
Implicit penultimate binary copolymerization model.
This model is a special case of the general (explicit) penultimate model, with a smaller number of independent parameters. As in the explicit version, the reactivity of a macroradical depends on the nature of the penultimate and terminal repeating units. A binary system is, thus, described by eight propagation reactions:
where \(k_{iii}\) are the homo-propagation rate coefficients and \(k_{ijk}\) are the cross-propagation coefficients. The six cross-propagation coefficients are specified via just four reactivity ratios, which are divided in two categories. There are two monomer reactivity ratios, which are defined as \(r_1=k_{111}/k_{112}=k_{211}/k_{212}\) and \(r_2=k_{222}/k_{221}=k_{122}/k_{121}\). Additionally, there are two radical reactivity ratios defined as \(s_1=k_{211}/k_{111}\) and \(s_2=k_{122}/k_{222}\). The latter influence the average propagation rate coefficient, but have no effect on the copolymer composition.
PARAMETER | DESCRIPTION |
---|---|
r1
|
Monomer reactivity ratio, \(r_1=k_{111}/k_{112}=k_{211}/k_{212}\).
TYPE:
|
r2
|
Monomer reactivity ratio, \(r_2=k_{222}/k_{221}=k_{122}/k_{121}\).
TYPE:
|
s1
|
Radical reactivity ratio, \(s_1=k_{211}/k_{111}\).
TYPE:
|
s2
|
Radical reactivity ratio, \(s_2=k_{122}/k_{222}\).
TYPE:
|
k1
|
Homopropagation rate coefficient of M1, \(k_1 \equiv k_{111}\).
TYPE:
|
k2
|
Homopropagation rate coefficient of M2, \(k_2 \equiv k_{222}\).
TYPE:
|
M1
|
Name of M1.
TYPE:
|
M2
|
Name of M2.
TYPE:
|
name
|
Name.
TYPE:
|
Source code in src/polykin/copolymerization/penultimate.py
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 |
|
F1 ¤
F1(
f1: Union[float, FloatArrayLike]
) -> Union[float, FloatArray]
Calculate the instantaneous copolymer composition, \(F_1\).
The calculation is handled by
inst_copolymer_binary
.
PARAMETER | DESCRIPTION |
---|---|
f1
|
Molar fraction of M1.
TYPE:
|
RETURNS | DESCRIPTION |
---|---|
float | FloatArray
|
Instantaneous copolymer composition, \(F_1\). |
Source code in src/polykin/copolymerization/terminal.py
85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 |
|
add_data ¤
add_data(
data: Union[CopoDataset, list[CopoDataset]]
) -> None
Add a copolymerization dataset for subsequent analysis.
PARAMETER | DESCRIPTION |
---|---|
data
|
Experimental dataset(s).
TYPE:
|
Source code in src/polykin/copolymerization/terminal.py
424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 |
|
azeotrope
property
¤
azeotrope: Optional[float]
Calculate the azeotrope composition.
An azeotrope (i.e., a point where \(F_1=f_1\)) only exists if both reactivity ratios are smaller than unity. In that case, the azeotrope composition is given by:
where \(r_1\) and \(r_2\) are the reactivity ratios.
RETURNS | DESCRIPTION |
---|---|
float | None
|
If an azeotrope exists, it returns its composition in terms of \(f_1\). |
drift ¤
drift(
f10: Union[float, FloatVectorLike],
x: Union[float, FloatVectorLike],
) -> FloatArray
Calculate drift of comonomer composition in a closed system for a given total monomer conversion.
In a closed binary system, the drift in monomer composition is given by the solution of the following differential equation:
with initial condition \(f_1(0)=f_{1,0}\), where \(f_1\) and \(F_1\) are, respectively, the instantaneous comonomer and copolymer composition of M1, and \(x\) is the total molar monomer conversion.
PARAMETER | DESCRIPTION |
---|---|
f10
|
Initial molar fraction of M1, \(f_{1,0}=f_1(0)\).
TYPE:
|
x
|
Value(s) of total monomer conversion values where the drift is to be evaluated.
TYPE:
|
RETURNS | DESCRIPTION |
---|---|
FloatArray
|
Monomer fraction of M1 at a given conversion, \(f_1(x)\). |
Source code in src/polykin/copolymerization/terminal.py
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 |
|
from_Qe
classmethod
¤
from_Qe(
Qe1: tuple[float, float],
Qe2: tuple[float, float],
k1: Optional[Arrhenius] = None,
k2: Optional[Arrhenius] = None,
M1: str = "M1",
M2: str = "M2",
name: str = "",
) -> TerminalModel
Construct TerminalModel
from Q-e values.
Alternative constructor that takes the \(Q\)-\(e\) values of the monomers as primary input instead of the reactivity ratios.
The conversion from Q-e to r is handled by
convert_Qe_to_r
.
PARAMETER | DESCRIPTION |
---|---|
Qe1
|
Q-e values of M1.
TYPE:
|
Qe2
|
Q-e values of M2.
TYPE:
|
k1
|
Homopropagation rate coefficient of M1, \(k_1 \equiv k_{11}\).
TYPE:
|
k2
|
Homopropagation rate coefficient of M2, \(k_2 \equiv k_{22}\).
TYPE:
|
M1
|
Name of M1.
TYPE:
|
M2
|
Name of M2.
TYPE:
|
name
|
Name.
TYPE:
|
Source code in src/polykin/copolymerization/terminal.py
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 |
|
kii ¤
kii(
f1: Union[float, FloatArray],
T: float,
Tunit: Literal["C", "K"] = "K",
) -> tuple[
Union[float, FloatArray], Union[float, FloatArray]
]
Pseudo-homopropagation rate coefficients.
In the implicit penultimate model, the pseudohomopropagation rate coefficients depend on the instantaneous comonomer composition according to:
where \(r_i\) are the monomer reactivity ratios, \(s_i\) are the radical reactivity ratios, and \(k_{iii}\) are the homopropagation rate coefficients.
PARAMETER | DESCRIPTION |
---|---|
f1
|
Molar fraction of M1.
TYPE:
|
T
|
Temperature. Unit =
TYPE:
|
Tunit
|
Temperature unit.
TYPE:
|
RETURNS | DESCRIPTION |
---|---|
tuple[float | FloatArray, float | FloatArray]
|
Tuple of pseudohomopropagation rate coefficients, (\(\bar{k}_{11}\), \(\bar{k}_{22}\)). |
Source code in src/polykin/copolymerization/penultimate.py
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 |
|
kp ¤
kp(
f1: Union[float, FloatArrayLike],
T: float,
Tunit: Literal["C", "K"] = "K",
) -> Union[float, FloatArray]
Calculate the average propagation rate coefficient, \(\bar{k}_p\).
The calculation is handled by
kp_average_binary
.
Note
This feature requires the attributes k11
and k22
to be defined.
PARAMETER | DESCRIPTION |
---|---|
f1
|
Molar fraction of M1.
TYPE:
|
T
|
Temperature. Unit =
TYPE:
|
Tunit
|
Temperature unit.
TYPE:
|
RETURNS | DESCRIPTION |
---|---|
float | FloatArray
|
Average propagation rate coefficient. Unit = L/(mol·s) |
Source code in src/polykin/copolymerization/terminal.py
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 |
|
plot ¤
plot(
kind: Literal["drift", "kp", "Mayo", "triads"],
show: Literal["auto", "all", "data", "model"] = "auto",
M: Literal[1, 2] = 1,
f0: Optional[Union[float, FloatVectorLike]] = None,
T: Optional[float] = None,
Tunit: Literal["C", "K"] = "K",
title: Optional[str] = None,
axes: Optional[Axes] = None,
return_objects: bool = False,
) -> Optional[tuple[Optional[Figure], Axes]]
Generate a plot of instantaneous copolymer composition, monomer composition drift, or average propagation rate coefficient.
PARAMETER | DESCRIPTION |
---|---|
kind
|
Kind of plot to be generated.
TYPE:
|
show
|
What informatation is to be plotted.
TYPE:
|
M
|
Index of the monomer to be used in input argument
TYPE:
|
f0
|
Initial monomer composition, \(f_i(0)\), as required for a monomer composition drift plot.
TYPE:
|
T
|
Temperature. Unit =
TYPE:
|
Tunit
|
Temperature unit.
TYPE:
|
title
|
Title of plot. If
TYPE:
|
axes
|
Matplotlib Axes object.
TYPE:
|
return_objects
|
If
TYPE:
|
RETURNS | DESCRIPTION |
---|---|
tuple[Figure | None, Axes] | None
|
Figure and Axes objects if return_objects is |
Source code in src/polykin/copolymerization/terminal.py
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 |
|
ri ¤
ri(
_,
) -> tuple[
Union[float, FloatArray], Union[float, FloatArray]
]
Return the evaluated reactivity ratios at the given conditions.
Source code in src/polykin/copolymerization/terminal.py
585 586 587 588 |
|
sequence ¤
sequence(
f1: Union[float, FloatArrayLike],
k: Optional[Union[int, IntArrayLike]] = None,
) -> dict[str, Union[float, FloatArray]]
Calculate the instantaneous sequence length probability or the number-average sequence length.
For a binary system, the probability of finding \(k\) consecutive units of monomer \(i\) in a chain is:
and the corresponding number-average sequence length is:
where \(P_{ii}\) is the transition probability \(i \rightarrow i\), which is a function of the monomer composition and the reactivity ratios.
References
- NA Dotson, R Galván, RL Laurence, and M Tirrel. Polymerization process modeling, Wiley, 1996, p. 177.
PARAMETER | DESCRIPTION |
---|---|
f1
|
Molar fraction of M1.
TYPE:
|
k
|
Sequence length, i.e., number of consecutive units in a chain.
If
TYPE:
|
RETURNS | DESCRIPTION |
---|---|
dict[str, float | FloatArray]
|
If |
Source code in src/polykin/copolymerization/terminal.py
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 |
|
transitions ¤
transitions(
f1: Union[float, FloatArrayLike]
) -> dict[str, Union[float, FloatArray]]
Calculate the instantaneous transition probabilities.
For a binary system, the transition probabilities are given by:
where \(i,j=1,2, i \neq j\).
References
- NA Dotson, R Galván, RL Laurence, and M Tirrel. Polymerization process modeling, Wiley, 1996, p. 178.
PARAMETER | DESCRIPTION |
---|---|
f1
|
Molar fraction of M1.
TYPE:
|
RETURNS | DESCRIPTION |
---|---|
dict[str, float | FloatArray]
|
Transition probabilities, {'11': \(P_{11}\), '12': \(P_{12}\), ... }. |
Source code in src/polykin/copolymerization/terminal.py
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 |
|
triads ¤
triads(
f1: Union[float, FloatArrayLike]
) -> dict[str, Union[float, FloatArray]]
Calculate the instantaneous triad fractions.
For a binary system, the triad fractions are given by:
where \(P_{ij}\) is the transition probability \(i \rightarrow j\) and \(F_i\) is the instantaneous copolymer composition.
References
- NA Dotson, R Galván, RL Laurence, and M Tirrel. Polymerization process modeling, Wiley, 1996, p. 179.
PARAMETER | DESCRIPTION |
---|---|
f1
|
Molar fraction of M1.
TYPE:
|
RETURNS | DESCRIPTION |
---|---|
dict[str, float | FloatArray]
|
Triad fractions, {'111': \(A_{111}\), '112': \(A_{112}\), '212': \(A_{212}\), ... }. |
Source code in src/polykin/copolymerization/terminal.py
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 |
|