yaeos__equilibria_boundaries_phase_envelopes_mp Module

Multiphase PT envelope calculation module.

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



Derived Types

type, public ::  PTEnvelMP

Multiphase PT envelope.

Components

Type Visibility Attributes Name Initial
type(MPPoint), public, allocatable :: points(:)

Array of converged points.

Type-Bound Procedures

procedure, public, nopass :: get_values_from_X
procedure, public, nopass :: solve_point
procedure, public :: write => write_envelope_PT_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.

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 pt_envelope(model, z, np, x_l0, w0, betas0, P0, T0, ns0, dS0, beta_w, points, max_pressure)

Calculation of a multiphase PT envelope.

Read more…

Arguments

Type IntentOptional Attributes Name
class(ArModel), intent(in) :: model
real(kind=pr), intent(in) :: z(:)

Mixture global composition.

integer, intent(in) :: np

Number of main phases.

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

Initial guess for the mole fractions of each phase. arranged as an array of size (np, nc), where nc is the number of components and np the number of main phases. Each row correspond to the composition of each main phaase.

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

Initial guess for the mole fractions of the reference/incipient phase.

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

Initial guess for the fractions of the main phases. arranged as an array of size (np), where np is the number of main phases.

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

Initial guess for the pressure [bar].

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

Initial guess for the temperature [K].

integer, intent(in) :: ns0

Number of the specified variable. The variable to be specified. This is the variable that will be used to calculate the first point of the envelope. The variable can be any of the variables in the vector X, but it is recommended to use the temperature or pressure. The variables are aranged as follows:

Read more…
real(kind=pr), intent(in) :: dS0

Step size of the specification for the next point. This is the step size that will be used to calculate the next point. Inside the algorithm this value is modified to adapt the step size to facilitate the convergence of each point.

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

Fraction of the reference (incipient) phase.

integer, intent(in), optional :: points

Number of points to calculate.

real(kind=pr), intent(in), optional :: max_pressure

Maximum pressure [bar] to calculate. If the pressure of the point is greater than this value, the calculation is stopped. This is useful to avoid calculating envelopes that go to infinite values of pressure.

Return Value type(PTEnvelMP)


Subroutines

public subroutine pt_F_NP(model, z, np, beta_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

Model to use.

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

Mixture global composition.

integer, intent(in) :: np

Number of main phases.

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

Fraction of the reference (incipient) 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 detect_critical(nc, np, X, dXdS, ns, dS, S)

Detect if the system is close to a critical point.

Read more…

Arguments

Type IntentOptional Attributes Name
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) :: 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.

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

Specification value.

private subroutine get_values_from_X(X, np, z, x_l, w, betas, P, T)

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) :: z(:)

Mixture composition.

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

Mole fractions of the main phases.

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

Mole fractions of the incipient phase.

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

Fractions of the main phases.

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

Pressure [bar].

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

Temperature [K].

private subroutine solve_point(model, z, np, beta_w, X, ns, S, dXdS, F, df, iters, max_iterations)

Arguments

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

Model to use.

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

Mixture global composition.

integer, intent(in) :: np

Number of main phases

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

Fraction of the reference (incipient) 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 to solve the current 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_PT_MP(env, unit)

Arguments

Type IntentOptional Attributes Name
class(PTEnvelMP), intent(in) :: env
integer, intent(in) :: unit