Attractive parameter mixing rule with quadratic mix.
Takes the all the pure components attractive parameters and their derivatives with respect to temperature and mix them with the Van der Waals quadratic mixing rule:
Inside the routine the matrix is calculated using the
procedure contained in the QMR
object, this procedures defaults
to the common combining rule:
The procedure can be overloaded by a common one that respects the interface get_aij
type(QMR) :: my_mixing_rule
my_mixing_rule%aij => new_aij_procedure
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(QMR), | intent(in) | :: | self |
Mixing rule object. |
||
real(kind=pr), | intent(in) | :: | n(:) |
Moles vector [mol] |
||
real(kind=pr), | intent(in) | :: | T |
Temperature [K] |
||
real(kind=pr), | intent(in) | :: | ai(:) |
Pure components attractive parameters |
||
real(kind=pr), | intent(in) | :: | daidt(:) |
|
||
real(kind=pr), | intent(in) | :: | daidt2(:) |
|
||
real(kind=pr), | intent(out) | :: | D |
Mixture attractive parameter |
||
real(kind=pr), | intent(out) | :: | dDdT |
|
||
real(kind=pr), | intent(out) | :: | dDdT2 |
|
||
real(kind=pr), | intent(out) | :: | dDi(:) |
|
||
real(kind=pr), | intent(out) | :: | dDidT(:) |
|
||
real(kind=pr), | intent(out) | :: | dDij(:,:) |
|
Type | Visibility | Attributes | Name | Initial | |||
---|---|---|---|---|---|---|---|
real(kind=pr), | public | :: | aij(size(ai),size(ai)) | ||||
real(kind=pr), | public | :: | aux | ||||
real(kind=pr), | public | :: | aux2 | ||||
real(kind=pr), | public | :: | daijdt(size(ai),size(ai)) | ||||
real(kind=pr), | public | :: | daijdt2(size(ai),size(ai)) | ||||
integer, | public | :: | i | ||||
integer, | public | :: | j | ||||
integer, | public | :: | nc |