yaeos__equilibria_boundaries_phase_envelopes_mp_tx Module

Multiphase Px envelope calculation module.

This module contains the functions to calculate the PT envelope of a mixture with multiple phases.



Derived Types

type, public ::  TXEnvelMP

Multiphase T envelope.

Components

Type Visibility Attributes Name Initial
real(kind=pr), public, allocatable :: Tc(:)

Critical temperatures [K]

real(kind=pr), public, allocatable :: ac(:)

Critical alphas.

real(kind=pr), public, allocatable :: alpha(:)

Array of values for each point in the envelope.

character(len=14), public, allocatable :: kind_w(:)

Kind of the reference phase

character(len=14), public, allocatable :: kinds_x(:)

Kinds of the main phases

type(MPPoint), public, allocatable :: points(:)

Array of converged points.

real(kind=pr), public, allocatable :: z0(:)

Initial mixture composition.

real(kind=pr), public, allocatable :: zi(:)

Second mixture composition.

Type-Bound Procedures

procedure, public, nopass :: get_values_from_X
procedure, public, nopass :: solve_point
procedure, public :: write => write_envelope_TX_MP

type, private ::  MPPoint

Multiphase equilibria point.

Components

Type Visibility Attributes Name Initial
real(kind=pr), public :: P

Pressure [bar]

real(kind=pr), public :: T

Temperature [K]

real(kind=pr), public :: beta_w

Fraction of the reference (incipient) phase.

real(kind=pr), public, allocatable :: betas(:)

Fractions of the main phases.

integer, public :: iters

Number of iterations needed to converge the point.

character(len=14), public :: kind_w

Kind of the reference phase.

character(len=14), public, allocatable :: kinds_x(:)

Kinds of the main phases.

integer, public :: nc

Number of components

integer, public :: np

Number of phases

integer, public :: ns

Number of the specified variable.

real(kind=pr), public, allocatable :: w(:)

Mole fractions of the incipient phase.

real(kind=pr), public, allocatable :: x_l(:,:)

Mole fractions of the main phases.


Functions

public function tx_envelope(model, z0, zi, np, P, kinds_x, kind_w, x_l0, w0, betas0, T0, alpha0, ns0, dS0, beta_w, points)

Calculate the multiphase isobar of a mixture.

Read more…

Arguments

Type IntentOptional Attributes Name
class(ArModel), intent(in) :: model

Model to use for the calculations.

real(kind=pr), intent(in) :: z0(:)

Initial mixture composition.

real(kind=pr), intent(in) :: zi(:)

Second mixture composition.

integer, intent(in) :: np

Number of main phases.

real(kind=pr), intent(in) :: P

Pressure [bar].

character(len=14), intent(in) :: kinds_x(np)

Kinds of the main phases. Can be “liquid”, “vapor”, “stable”

character(len=14), intent(in) :: kind_w

Kind of the reference phase. Can be “liquid”, “vapor”, “stable”

real(kind=pr), intent(in) :: x_l0(np,size(z0))

Mole fractions matrix of the main phases at the initial point. Each row corresponds to a phase, and each column corresponds to a component.

real(kind=pr), intent(in) :: w0(size(z0))

Mole fractions of the incipient/reference phase at the initial point.

real(kind=pr), intent(in) :: betas0(np)

Mole fractions of each main phases at the initial point.

real(kind=pr), intent(in) :: T0

Initial Temperature [K].

real(kind=pr), intent(in) :: alpha0

Initial value of the variable.

integer, intent(in) :: ns0

Number of the specified variable. This is the variable that will be used to specify the next point. the first variables corresponds to the of each phase, sorted by phases ([\ln K_1^1, \ln K_2^1, \dots, \ln K_{nc}^1, \dots, \ln K_{nc}^{np}]). The next variables are the mole fractions of each main phase. The last two variables are the temperature and the variable. Respectively

real(kind=pr), intent(in) :: dS0

Step size of the specification for the next point.

real(kind=pr), intent(in) :: beta_w

Fraction of the reference (incipient) phase.

integer, intent(in), optional :: points

Number of points to calculate in the envelope. If not specified, the default value is 1000.

Return Value type(TXEnvelMP)


Subroutines

public subroutine tx_F_NP(model, z0, zi, np, P, beta_w, kinds_x, kind_w, X, ns, S, F, df)

Function to solve at each point of a multi-phase envelope.

Arguments

Type IntentOptional Attributes Name
class(ArModel), intent(in) :: model
real(kind=pr), intent(in) :: z0(:)
real(kind=pr), intent(in) :: zi(:)
integer, intent(in) :: np

Number of main phases.

real(kind=pr), intent(in) :: P

Pressure [bar]

real(kind=pr), intent(in) :: beta_w

Fraction of the reference (incipient) phase.

character(len=14), intent(in) :: kinds_x(np)

Kinds of the main phases.

character(len=14), intent(in) :: kind_w

Kind of the reference phase.

real(kind=pr), intent(in) :: X(:)

Vector of variables.

integer, intent(in) :: ns

Number of specification.

real(kind=pr), intent(in) :: S

Specification value.

real(kind=pr), intent(out) :: F(size(X))

Vector of functions valuated.

real(kind=pr), intent(out) :: df(size(X),size(X))

Jacobian matrix.

private subroutine get_values_from_X(X, np, z0, zi, beta_w, x_l, w, betas, T, alpha)

Extract the values of the variables from the vector X.

Read more…

Arguments

Type IntentOptional Attributes Name
real(kind=pr), intent(in) :: X(:)

Vector of variables.

integer, intent(in) :: np

Number of main phases.

real(kind=pr), intent(in) :: z0(:)

Initial mixture composition.

real(kind=pr), intent(in) :: zi(:)

Second mixture composition.

real(kind=pr), intent(in) :: beta_w

Reference phase beta.

real(kind=pr), intent(out) :: x_l(np,size(z0))

Mole fractions of the main phases.

real(kind=pr), intent(out) :: w(size(z0))

Mole fractions of the incipient phase.

real(kind=pr), intent(out) :: betas(np)

Fractions of the main phases.

real(kind=pr), intent(out) :: T

Pressure [bar].

real(kind=pr), intent(out) :: alpha

.

private subroutine solve_point(model, z0, zi, np, P, beta_w, kinds_x, kind_w, X, ns, S, dXdS, F, df, iters, max_iterations)

Solve a point in the multiphase envelope.

Read more…

Arguments

Type IntentOptional Attributes Name
class(ArModel), intent(in) :: model

Model to use for the calculations.

real(kind=pr), intent(in) :: z0(:)

Initial mixture composition.

real(kind=pr), intent(in) :: zi(:)

Second mixture composition.

integer, intent(in) :: np

Number of main phases.

real(kind=pr), intent(in) :: P

Presure [bar].

real(kind=pr), intent(in) :: beta_w

Fraction of the reference (incipient) phase

character(len=14), intent(in) :: kinds_x(np)

Kinds of the main phases

character(len=14), intent(in) :: kind_w

Kind of the reference phase

real(kind=pr), intent(inout) :: X(:)

Vector of variables

integer, intent(in) :: ns

Number of specification

real(kind=pr), intent(in) :: S

Specification value

real(kind=pr), intent(in) :: dXdS(size(X))
real(kind=pr), intent(out) :: F(size(X))

Vector of functions valuated

real(kind=pr), intent(out) :: df(size(X),size(X))

Jacobian matrix

integer, intent(out) :: iters

Number of iterations needed to converge the point.

integer, intent(in) :: max_iterations

Maximum number of iterations to solve the point.

private subroutine update_specification(its, nc, np, X, dF, dXdS, ns, dS)

Change the specified variable for the next step.

Read more…

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: its

Iterations to solve the current point.

integer, intent(in) :: nc

Number of components in the mixture.

integer, intent(in) :: np

Number of main phases.

real(kind=pr), intent(inout) :: X(:)

Vector of variables.

real(kind=pr), intent(inout) :: dF(:,:)

Jacobian matrix.

real(kind=pr), intent(inout) :: dXdS(:)

Sensitivity of the variables wrt the specification.

integer, intent(inout) :: ns

Number of the specified variable.

real(kind=pr), intent(inout) :: dS

Step size of the specification for the next point.

private subroutine write_envelope_TX_MP(env, unit)

Write the multiphase envelope to a file.

Read more…

Arguments

Type IntentOptional Attributes Name
class(TXEnvelMP), intent(in) :: env

Envelope to write.

integer, intent(in) :: unit

Unit to write the envelope to.