yaeos__equilibria_boundaries_phase_envelopes_pt3 Module



Variables

Type Visibility Attributes Name Initial
real(kind=pr), public, parameter :: lnK_min = 2.0_pr

Derived Types

type, public ::  PTEnvel3

Components

Type Visibility Attributes Name Initial
real(kind=pr), public, allocatable :: P(:)
real(kind=pr), public, allocatable :: S(:)
real(kind=pr), public, allocatable :: T(:)
real(kind=pr), public, allocatable :: beta(:)

Mole fraction between phase x and phase y

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

Mole fraction of phase x

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

Mole fraction of phase x

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

Mole fraction of phase x


Functions

public function critical_interpol(Xnew, Xold, idx) result(a)

Critical point interpolation

Read more…

Arguments

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

New value of the variables

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

Old value of the variables

integer, intent(in) :: idx(:)

Index of the variables to interpolate

Return Value real(kind=pr)

public function pt_envelope_3ph(model, z, x0, y0, w0, beta0, P0, T0, ns0, dS0, points) result(envelope)

Arguments

Type IntentOptional Attributes Name
class(ArModel), intent(in) :: model
real(kind=pr), intent(in) :: z(:)
real(kind=pr), intent(in) :: x0(:)
real(kind=pr), intent(in) :: y0(:)
real(kind=pr), intent(in) :: w0(:)
real(kind=pr), intent(in) :: beta0
real(kind=pr), intent(in) :: P0
real(kind=pr), intent(in) :: T0
integer, intent(in) :: ns0
real(kind=pr), intent(in) :: dS0
integer, intent(in) :: points

Return Value type(PTEnvel3)


Subroutines

public subroutine detect_critical(X, dXdS, ns, S, dS)

Critical point detection

Read more…

Arguments

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

Vector of variables

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

Variation of variables wrt S

integer, intent(inout) :: ns

Number of specified variable

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

Specification value

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

Step in specification

public subroutine get_values_from_X(z, Xvars, x, y, w, P, T, beta)

Arguments

Type IntentOptional Attributes Name
real(kind=pr), intent(in) :: z(:)
real(kind=pr), intent(in) :: Xvars(size(z)*2+3)
real(kind=pr), intent(out) :: x(size(z))
real(kind=pr), intent(out) :: y(size(z))
real(kind=pr), intent(out) :: w(size(z))
real(kind=pr), intent(out) :: P
real(kind=pr), intent(out) :: T
real(kind=pr), intent(out) :: beta

public subroutine pt_F_three_phases(model, z, Xvars, ns, S, F, df)

Function to solve at each point of a three phase envelope.

Read more…

Arguments

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

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(Xvars))

Vector of functions valuated

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

Jacobian matrix

public subroutine solve_point(model, z, ns, S, X, F, dF, its, maxits)

Arguments

Type IntentOptional Attributes Name
class(ArModel), intent(in) :: model
real(kind=pr), intent(in) :: z(:)
integer, intent(in) :: ns
real(kind=pr), intent(in) :: S
real(kind=pr), intent(inout) :: X(:)
real(kind=pr), intent(out) :: F(:)
real(kind=pr), intent(out) :: dF(:,:)
integer, intent(inout) :: its
integer, intent(in) :: maxits

public subroutine update_specification(its, X, dF, dXdS, ns, dS)

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: its
real(kind=pr), intent(inout) :: X(:)
real(kind=pr), intent(inout) :: dF(:,:)
real(kind=pr), intent(inout) :: dXdS(:)
integer, intent(inout) :: ns
real(kind=pr), intent(inout) :: dS