pt_envelope Function

public function pt_envelope(model, z, np, kinds_x, kind_w, x_l0, w0, betas0, P0, T0, ns0, dS0, beta_w, points, max_pressure)

Uses

  • proc~~pt_envelope~~UsesGraph proc~pt_envelope pt_envelope module~yaeos__auxiliar yaeos__auxiliar proc~pt_envelope->module~yaeos__auxiliar module~yaeos__constants yaeos__constants module~yaeos__auxiliar->module~yaeos__constants iso_fortran_env iso_fortran_env module~yaeos__constants->iso_fortran_env

pt_envelope

Calculation of a multiphase PT envelope.

Description

Calculates a PT envelope is calculated using the continuation method. The envelope is calculated by solving the system of equations for each point of the envelope. The system of equations is solved using the Newton-Raphson method.

This function requires the system specification conditions, which are the fluid composition (\z), the number of phases that are not incipient; defined as , proper intialization values, the variables that end with 0 are the initial guess; the mole fraction of the reference phase beta_w which when it is equal to 0 means that we are calculating a phase boundary.

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.

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

Kind of the main phases.

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

Kind of the reference phase.

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:

  • X(1:nc*np) = ln(K_i^l):
  • X(nc*np+1:nc*np+np) = \beta_i^l: Fraction of each main phase.
  • X(nc*np+np+1) = ln(P): Pressure [bar].
  • X(nc*np+np+2) = ln(T): Temperature [K].
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)


Calls

proc~~pt_envelope~~CallsGraph proc~pt_envelope pt_envelope interface~optval optval proc~pt_envelope->interface~optval proc~check_critical_jump check_critical_jump proc~pt_envelope->proc~check_critical_jump proc~detect_critical~3 detect_critical proc~pt_envelope->proc~detect_critical~3 proc~get_values_from_x~3 PTEnvelMP%get_values_from_X proc~pt_envelope->proc~get_values_from_x~3 proc~solve_point~4 PTEnvelMP%solve_point proc~pt_envelope->proc~solve_point~4 proc~update_specification~3 update_specification proc~pt_envelope->proc~update_specification~3 proc~optval_character optval_character interface~optval->proc~optval_character proc~optval_integer optval_integer interface~optval->proc~optval_integer proc~optval_real optval_real interface~optval->proc~optval_real proc~interpol interpol proc~check_critical_jump->proc~interpol proc~detect_critical~3->proc~interpol eye eye proc~solve_point~4->eye proc~pt_f_np pt_F_NP proc~solve_point~4->proc~pt_f_np proc~solve_system solve_system proc~solve_point~4->proc~solve_system proc~update_specification~3->proc~solve_system dvnvdb dvnvdb proc~pt_f_np->dvnvdb dvnvdlnkl dvnvdlnkl proc~pt_f_np->dvnvdlnkl dvnvdn dvnvdn proc~pt_f_np->dvnvdn dvwdb dvwdb proc~pt_f_np->dvwdb dvwdlnkl dvwdlnkl proc~pt_f_np->dvwdlnkl proc~lnphi_pt ArModel%lnphi_pt proc~pt_f_np->proc~lnphi_pt none~dgesv dgesv proc~solve_system->none~dgesv proc~lnphi_vt ArModel%lnphi_vt proc~lnphi_pt->proc~lnphi_vt proc~volume~3 ArModel%volume proc~lnphi_pt->proc~volume~3 residual_helmholtz residual_helmholtz proc~lnphi_vt->residual_helmholtz get_v0 get_v0 proc~volume~3->get_v0 interface~newton newton proc~volume~3->interface~newton proc~newton_1d newton_1d interface~newton->proc~newton_1d

Variables

Type Visibility Attributes Name Initial
real(kind=pr), private :: F(size(z)*np+np+2)

Vector of functions valuated.

real(kind=pr), private :: P

Pressure [bar].

real(kind=pr), private :: Pc

Critical pressure [bar]

real(kind=pr), private :: S

Specified value

real(kind=pr), private :: T

Temperature [K].

real(kind=pr), private :: Tc

Critical temperature [K]

real(kind=pr), private :: Vl(np)

Molar volumes of the main phases [L/mol]

real(kind=pr), private :: Vw

Molar volume of the reference phase [L/mol]

real(kind=pr), private :: X(size(z)*np+np+2)

Vector of variables.

real(kind=pr), private :: X0(size(X))

Initial guess for the point

real(kind=pr), private :: X_last_converged(size(X))

Last converged point

real(kind=pr), private :: Xc(size(X))

Vector of variables at the critical point

real(kind=pr), private :: betas(np)

Fractions of the main phases.

real(kind=pr), private :: dF(size(z)*np+np+2,size(z)*np+np+2)

Jacobian matrix.

real(kind=pr), private :: dS

Step size of the specification for the next point

real(kind=pr), private :: dX(size(z)*np+np+2)

Step for next point estimation.

real(kind=pr), private :: dXdS(size(z)*np+np+2)

Sensitivity of the variables wrt the specification.

type(MPPoint), private, allocatable :: env_points(:)

Array of converged points.

logical, private :: found_critical

If true, a critical point was found

integer, private :: i

Point calculation index

integer, private :: iP

Index of the pressure variable.

integer, private :: iT

Index of the temperature variable.

integer, private :: inner

Number of times a failed point is retried to converge

integer, private :: its

Number of iterations to solve the current point.

logical, private :: jumped_critical

If true, a critical point was jumped

integer, private :: l

Phase index

integer, private :: lb

Lower bound, index of the first component of a phase

real(kind=pr), private :: max_P

Maximum pressure [bar] to calculate.

integer, private :: max_iterations = 100

Maximum number of iterations to solve the point.

integer, private :: nc

Number of components.

integer, private :: ns

Number of the specified variable

integer, private :: number_of_points

Number of points to calculate.

type(MPPoint), private :: point

Converged point.

integer, private :: ub

Upper bound, index of the last component of a phase

real(kind=pr), private :: w(size(z))

Mole fractions of the incipient phase.

character(len=14), private :: w_kind
character(len=14), private :: x_kinds(np)
real(kind=pr), private :: x_l(np,size(z))

Mole fractions of the main phases.