conmin_module Module

CONMIN is a Fortran package for the solution of linear or nonlinear constrained optimization problems. The basic optimization algorithm is the Method of Feasible Directions. The user must provide a main calling program and an external routine to evaluate the objective and constraint functions and to provide gradient information. If analytic gradients of the objective or constraint functions are not available, this information is calculated by finite difference. While the program is intended primarily for efficient solution of constrained problems, unconstrained function minimization problems may also be solved, and the conjugate direction method of Fletcher and Reeves is used for this purpose.

History

  • BY G. N. VANDERPLAATS NASA-AMES RESEARCH CENTER, MOFFETT FIELD, CALIF.
  • Modernized by Jacob Williams, March 2025.

See also

  • CONMIN - A FORTRAN PROGRAM FOR CONSTRAINED FUNCTION MINIMIZATION: USER'S MANUAL, BY G. N. VANDERPLAATS, NASA TM X-62,282, AUGUST, 1973.
  • conmin.zip from Alan Miller.

Uses

  • module~~conmin_module~~UsesGraph module~conmin_module conmin_module iso_fortran_env iso_fortran_env module~conmin_module->iso_fortran_env

Variables

Type Visibility Attributes Name Initial
integer, private, parameter :: wp = real64

Real working precision if not specified [8 bytes]

integer, public, parameter :: conmin_wp = wp

Working real precision


Abstract Interfaces

abstract interface

  • private subroutine report_f(me, iter, x, obj, g)

    Report function to call once per iteration to report the solution.

    Arguments

    Type IntentOptional Attributes Name
    class(conmin_class), intent(inout) :: me
    integer, intent(in) :: iter

    Iteration number

    real(kind=wp), intent(in), dimension(:) :: x

    Decision variables

    real(kind=wp), intent(in) :: obj

    Objective function value

    real(kind=wp), intent(in), dimension(:) :: g

    Constraint functions


Derived Types

type, public ::  conmin_class

main class.

Components

Type Visibility Attributes Name Initial
real(kind=wp), public :: small = 1.0e-20_wp

Small number

real(kind=wp), public :: large = 1.0e+20_wp

Large number

real(kind=wp), public :: delfun = 0.001_wp

Minimum relative change in the objective function to indicate convergence. If in ITRM consecutive iterations, ABS(1.0-OBJ(J-1)/OBJ(J))<DELFUN and the current design is feasible (all G(J)<=ABS(CT)), the minimization process is terminated. If the current design is infeasible (some G(J)>ABS(CT)), five iterations are required to terminate and this situation indicates that a feasible design may not exist.

real(kind=wp), public :: dabfun = 0.0_wp

Default value = 0.001 times the initial function value. Same as DELFUN except comparison is on absolute change in the objective function, ABS(OBJ(J)-OBJ(J-1)), instead of relative change.

real(kind=wp), public :: fdch = 0.01_wp

Not used if NFDG = 0. Relative change in decision variable X(I) in calculating finite difference gradients. For example, FDCH = 0.01 corresponds to a finite difference step of one percent of the value of the decision variable.

real(kind=wp), public :: fdchm = 0.01_wp

Not used if NFDG = 0. Minimum absolute step in finite difference gradient calculations. FDCHM applies to the unscaled variable values.

real(kind=wp), public :: ct = -0.1_wp

Not used if NCON = NSIDE = 0. Constraint thickness parameter. If CT<=G(J)<=ABS(CT), G(J) is defined as active. If G(J)>ABS(CT), G(J) is said to be violated. If G(J)<CT, G(J) is not active. CT is sequentially reduced in magnitude during the optimization process. If ABS(CT) is very small, one or more constraints may be active on one iteration and inactive on the next, only to become active again on a subsequent iteration. This is often referred to as "zigzagging" between constraints. A wide initial value of the constraint thickness is desirable for highly nonlinear problems so that when a constraint becomes active it tends to remain active, thus reducing the zigzagging problem. The default value is usually adequate.

real(kind=wp), public :: ctmin = 0.004_wp

Not used if NCON = NSIDE = 0. Minimum absolute value of CT considered in the optimization process. CTMIN may be considered as "numerical zero" since it may not be meaningful to compare numbers smaller than CTMIN. The value of CTMIN is chosen to indicate that satisfaction of a constraint within this tolerance is acceptable. The default value is usually adequate.

real(kind=wp), public :: ctl = -0.01_wp

Not used if NCON = NSIDE = 0. Constraint thickness parameter for linear and side constraints. CTL is smaller in magnitude than CT because the zigzagging problem is avoided with linear and side constraints. The default value is usually adequate.

real(kind=wp), public :: ctlmin = 0.001_wp

Not used if NCON = NSIDE = 0. Minimum absolute value of CTL considered in the optimization process. The default value is usually adequate.

real(kind=wp), public :: alphax = 0.1_wp

the maximum fractional change in any component of X as an initial estimate for ALPHA in the one-dimensional search. That is, the initial ALPHA will be such that no component of X is changed by more than this amount. This only applies to those X(i) of magnitude greater than 0.1. If an optimization run shows numerous ALPHA = 0 results for the one-dimensional search, it may help to try ALPHAX less than the default. ALPHAX is changed by CONMIN depending on the progress of the optimization.

real(kind=wp), public :: abobj1 = 0.1_wp

the fractional change attempted as a first step in the one-dimensional search and is based on a linear approximation. ABOBJ1 is updated during the optimization, depending on progress. The initial step in the one-dimensional search is taken as the amount necessary to change OBJ by ABOBJ1*ABS(OBJ) or to change some X(i) by ALPHAX*ABS( X(i) ), whichever is less.

real(kind=wp), public :: theta = 1.0_wp

Not used if NCON = NSIDE = 0. Mean value of the push-off factor in the method of feasible directions. A larger value of THETA is desirable if the constraints, G(J), are known to be highly nonlinear, and a smaller value may be used if all G(J) are known to be nearly linear. The actual value of the push-off factor used in the program is a quadratic function of each G(J), varying from 0.0 for G(J) = CT to 4.0*THETA for G(J) = ABS(CT). A value of THETA = 0.0 is used in the program for constraints which are identified by the user to be strictly linear. THETA is called a "push-off" factor because it pushes the design away from the active constraints into the feasible region. The default value is usually adequate.

real(kind=wp), public :: obj

Value of objective function for the current decision variables, X(I), I = 1, NDV contained in vector X. Calculate OBJ if INFO = 1 or INFO = 2.

integer, public :: ndv

Number of decision variables, X(I), contained in vector X.

integer, public :: ncon

Number of constraint functions, G(J). NCON may be zero.

integer, public :: nside

Side constraint parameter. NSIDE = 0 signifies that the variables X(I) do not have lower or upper bounds. NSIDE>0 signifies that all variables X(I) have lower and upper bounds defined by VLB(I) and VUB(I) respectively. If one or more variables are not bounded while others are, the values of the lower and upper bounds on the unbounded variables must be taken as very large negative and positive values respectively (i.e., VLB(I) = -1.0E+10, VUB(I) = 1.0E+10).

integer, public :: iunit = output_unit

Unit number for printing. Should be non-zero.

integer, public :: iprint = 0

Print control. All printing is done on unit number iunit.

Read more…
integer, public :: nfdg

the finite difference gradient parameter:

Read more…
integer, public :: nscal

Scaling control parameter. The decision variables will be scaled linearly.

Read more…
integer, public :: linobj = 0

Not used if NCON = NSIDE = 0. Linear objective function identifier. If the objective, OBJ, is specifically known to be a strictly linear function of the decision variables, X(I), set LINOBJ = 1. If OBJ is a general nonlinear function, set LINOBJ = 0.

integer, public :: itmax = 10

Maximum number of iterations in the minimization process. If NFDG==0 each iteration requires one set of gradient computations (INFO = 3 or 4) and approximately three function evaluations (INFO = 1 or 2). If NFDG>0 each iteration requires approximately NDV + 3 function evaluations (INFO = 1 or 2).

integer, public :: itrm = 3

Number of consecutive iterations to indicate convergence by relative or absolute changes, DELFUN or DABFUN.

integer, public :: icndir

Default value = NDV + 1. Conjugate direction restart parameter. If the function is currently unconstrained, (all G(J)<CT or NCON = NSIDE = 0), Fletcher-Reeves conjugate direction method will be restarted with a steepest descent direction every ICNDIR iterations. If ICNDIR = 1 only steepest descent will be used.

integer, public :: igoto = 0

Reverse communication flag. CONMIN executes according to the parameter IGOTO which must be initialized to zero.

integer, public :: nac

Number of active and violated constraints (G(J)>=CT). Calculate NAC if INFO = 4 and NFDG = 0.

integer, public :: info Read more…
integer, public :: infog Read more…
integer, public :: iter

Iteration number. The optimization process is iterative so that the vector of decision variables at the Kth iteration is defined by X(K) = X(K - 1) + ALPHA*S(K), where in this case K refers to the iteration number and the components X(I) are all changed simultaneously. ALPHA is defined as the move parameter and is printed if the print control IPRINT>=3. S is the move direction.

real(kind=wp), public :: dct
real(kind=wp), public :: dctl
real(kind=wp), public :: abobj
real(kind=wp), public :: cta
real(kind=wp), public :: ctam
real(kind=wp), public :: ctbm
real(kind=wp), public :: obj1
real(kind=wp), public :: slope
real(kind=wp), public :: dx
real(kind=wp), public :: dx1
real(kind=wp), public :: fi
real(kind=wp), public :: xi
real(kind=wp), public :: dftdf1
real(kind=wp), public :: alp
real(kind=wp), public :: fff
real(kind=wp), public :: a1
real(kind=wp), public :: a2
real(kind=wp), public :: a3
real(kind=wp), public :: a4
real(kind=wp), public :: f1
real(kind=wp), public :: f2
real(kind=wp), public :: f3
real(kind=wp), public :: f4
real(kind=wp), public :: cv1
real(kind=wp), public :: cv2
real(kind=wp), public :: cv3
real(kind=wp), public :: cv4
real(kind=wp), public :: app
real(kind=wp), public :: alpca
real(kind=wp), public :: alpfes
real(kind=wp), public :: alpln
real(kind=wp), public :: alpmin
real(kind=wp), public :: alpnc
real(kind=wp), public :: alpsav
real(kind=wp), public :: alpsid
real(kind=wp), public :: alptot
real(kind=wp), public :: rspace
real(kind=wp), public :: phi = 5.0_wp

Not used if NCON = NSIDE = 0. Participation coefficient, used if a design is infeasible (one or more G(J)>ABS(CT)). PHI is a measure of how hard the design will be "pushed" towards the feasible region and is, in effect, a penalty parameter. If in a given problem, a feasible solution cannot be obtained with the default value, PHI should be increased, and the problem run again. If a feasible solution cannot be obtained with PHI = 100, it is probable that no feasible solution exists. The default value is usually adequate.

integer, public :: jdir
integer, public :: iobj
integer, public :: kobj
integer, public :: kcount
integer, public :: nfeas
integer, public :: mscal
integer, public :: ncobj
integer, public :: nvc
integer, public :: kount
integer, public :: icount
integer, public :: igood1
integer, public :: igood2
integer, public :: igood3
integer, public :: igood4
integer, public :: ibest
integer, public :: iii
integer, public :: nlnc
integer, public :: jgoto
integer, public :: ispace(2)
integer, public :: ncal(2)

Bookkeeping information. NCAL(1) gives the number of times external routine SUB1 was called with INFO = 1. NCAL(2) gives the number of times INFO = 2. NCAL(3) gives the number of times INFO = 3 and NCAL(4) gives the number of times INFO = 4.

procedure(report_f), public, pointer :: report => null()

a function to call once per iteration to report the solution

Type-Bound Procedures

procedure, public :: solve => conmin
procedure, private :: cnmn01
procedure, private :: cnmn02
procedure, private :: cnmn03
procedure, private :: cnmn04
procedure, private :: cnmn05
procedure, private :: cnmn06
procedure, private :: cnmn07

Subroutines

private subroutine conmin(me, x, vlb, vub, g, scal, df, a, s, g1, g2, b, c, isc, ic, ms1, n1, n2, n3, n4, n5)

Routine to solve constrained or unconstrained function minimization.

Read more…

Arguments

Type IntentOptional Attributes Name
class(conmin_class), intent(inout) :: me
real(kind=wp), intent(inout) :: x(n1)

Vector of decision variables, X(I), I = 1, NDV. The initial X-vector contains the user's best estimate of the set of optimum design variables.

real(kind=wp), intent(inout) :: vlb(n1)

Used only if NSIDE/=0. VLB(I) is the lower allowable value (lower bound) of variable X(I). If one or more variables, X(I), do not have lower bounds, the corresponding VLB(I) must be initialized to a very large negative number (say -1.0E+10).

real(kind=wp), intent(inout) :: vub(n1)

Used only if NSIDE/=0. VUB(I) is the maximum allowable value (upper bound) of X(I). If one or more variables, X(I), do not have upper bounds, the corresponding VUB(I) must be initialized to a very large positive number (say 1.0E+10).

real(kind=wp), intent(inout) :: g(n2)

Not used if NCON = NSIDE = 0. Vector containing all constraint functions, G(J), J = 1, NCON for current decision variables, X. Calculate G(J), J = 1, NCON if INFO = 2.

real(kind=wp), intent(inout) :: scal(n1)

Not used if NSCAL = 0. Vector of scaling parameters. If NSCAL>0 vector SCAL need not be initialized since SCAL will be defined in CONMIN and its associated routines. If NSCAL<0, vector SCAL is initialized in the main program, and the scaled variables X(I) = X(I)/SCAL(I). Efficiency of the optimization process can sometimes be improved if the variables are either normalized or are scaled in such a way that the partial deri- vative of the objective function, OBJ, with respect to variable X(I) is of the same order of magnitude for all X(I). SCAL(I) must be greater than zero because a negative value of SCAL(I) will result in a change of sign of X(I) and possibly yield erroneous optimization results. The decision of if, and how, the variables should be scaled is highly problem dependent, and some experimentation is desirable for any given class of problems.

real(kind=wp), intent(inout) :: df(n1)

Analytic gradient of the objective function for the current decision variables, X(I). DF(I) contains the partial derivative of OBJ with respect to X(I). Calculate DF(I), I = 1, NDV if INFO = 3 or INFO = 4 and if NFDG = 0 or NFDG = 2.

real(kind=wp), intent(inout) :: a(n1,n3)

Not used if NCON = NSIDE = 0. Gradients of active or violated constraints, for current decision variables, X(I). A(J,I) contains the gradient of the Jth active or violated constraint, G(J), with respect to the Ith decision variable, X(I) for J = 1, NAC and I = 1, NDV. Calculate if INFO = 4 and NFDG = 0.

real(kind=wp), intent(inout) :: s(n1)

Move direction in the NDV-dimensional optimization space. S(I) gives the rate at which variable X(I) changes with respect to ALPHA.

real(kind=wp), intent(inout) :: g1(n2)

Not used if NCON = NSIDE = NSCAL = 0. Used for temporary storage of constraint values G(J), J = 1, NCON and decision variables X(I), I = 1, NDV.

real(kind=wp), intent(inout) :: g2(n2)

Not used if NCON = NSIDE = 0. Used for temporary storage of constraint values G(J), J = 1, NCON.

real(kind=wp), intent(inout) :: b(n3,n3)

Not used if NCON = NSIDE = 0. Used in determining direction vector S for constrained minimization problems. Array B may be used for temporary storage in external routine SUB1.

real(kind=wp), intent(inout) :: c(n4)

Not used in NCON = NSIDE = 0. Used with array B in determining direction vector S for constrained minimization problems. Used for temporary storage of vector X if NSCAL/=0. routine SUB1.

integer, intent(inout) :: isc(n2)

Not used if NCON = 0. Linear constraint identification vector. If constraint G(J) is known to be a linear function of the decision variables, X(I), ISC(I) should be initialized to ISC(I) = 1. If constraint G(J) is nonlinear ISC(I) is initialized to ISC(I) = 0. Identification of linear constraints may improve efficiency of the optimization process and is therefore desirable, but is not essential. If G(J) is not specifically known to be linear, set ISC(I) = 0.

integer, intent(inout) :: ic(n3)

Identifies which constraints are active or violated. IC(J) contains the number of the Jth active or violated constraint for J = 1, NAC. For example, if G(10) is the first active or violated constraint (G(J)<CT, J = 1,9), set IC(1) = 10. Calculate if INFO = 4 and NFDG = 0.

integer, intent(inout) :: ms1(n5)

Not used if NCON = NSIDE = 0. Used with array B in determining direction vector S for constrained minimization problems. Array MS1 may be used for temporary storage in external routine SUB1.

integer, intent(in) :: n1

N1 = NDV + 2

integer, intent(in) :: n2

N2 = NCON + 2*NDV

integer, intent(in) :: n3

N3 = NACMX1

integer, intent(in) :: n4

N4 = MAX (N3,NDV)

integer, intent(in) :: n5

N5 = 2*N4

private subroutine cnmn01(me, jgoto, x, df, g, isc, ic, a, g1, vub, scal, ncal, dx, dx1, fi, xi, iii, n1, n2, n3)

Routine to calculate gradient information by finite difference.

Read more…

Arguments

Type IntentOptional Attributes Name
class(conmin_class), intent(inout) :: me
integer, intent(inout) :: jgoto
real(kind=wp), intent(inout) :: x(:)
real(kind=wp), intent(inout) :: df(:)
real(kind=wp), intent(inout) :: g(n2)
integer, intent(in) :: isc(n2)
integer, intent(inout) :: ic(n3)
real(kind=wp), intent(inout) :: a(n1,n3)
real(kind=wp), intent(inout) :: g1(n2)
real(kind=wp), intent(in) :: vub(:)
real(kind=wp), intent(in) :: scal(:)
integer, intent(inout) :: ncal(2)
real(kind=wp), intent(out) :: dx
real(kind=wp), intent(out) :: dx1
real(kind=wp), intent(out) :: fi
real(kind=wp), intent(out) :: xi
integer, intent(out) :: iii
integer, intent(in) :: n1
integer, intent(in) :: n2
integer, intent(in) :: n3

private subroutine cnmn02(me, ncalc, slope, dftdf1, df, s)

Routine to determine conjugate direction vector or direction of steepest descent for unconstrained function minimization.

Read more…

Arguments

Type IntentOptional Attributes Name
class(conmin_class), intent(inout) :: me
integer, intent(inout) :: ncalc

NCALC = CALCULATION CONTROL.

Read more…
real(kind=wp), intent(out) :: slope
real(kind=wp), intent(inout) :: dftdf1
real(kind=wp), intent(in) :: df(:)
real(kind=wp), intent(inout) :: s(:)

private subroutine cnmn03(me, x, s, slope, alp, fff, a1, a2, a3, a4, f1, f2, f3, f4, app, ncal, kount, jgoto)

Routine to solve one-dimensional search in unconstrained minimization using 2-point quadratic interpolation, 3-point cubic interpolation and 4-point cubic interpolation.

Read more…

Arguments

Type IntentOptional Attributes Name
class(conmin_class), intent(inout) :: me
real(kind=wp), intent(inout) :: x(:)
real(kind=wp), intent(inout) :: s(:)
real(kind=wp), intent(inout) :: slope

SLOPE = INITIAL FUNCTION SLOPE = S-TRANSPOSE TIMES DF. SLOPE MUST BE NEGATIVE.

real(kind=wp), intent(inout) :: alp

PROPOSED MOVE PARAMETER

real(kind=wp), intent(inout) :: fff
real(kind=wp), intent(inout) :: a1
real(kind=wp), intent(inout) :: a2
real(kind=wp), intent(inout) :: a3
real(kind=wp), intent(inout) :: a4
real(kind=wp), intent(inout) :: f1
real(kind=wp), intent(inout) :: f2
real(kind=wp), intent(inout) :: f3
real(kind=wp), intent(inout) :: f4
real(kind=wp), intent(inout) :: app
integer, intent(inout) :: ncal(2)
integer, intent(out) :: kount
integer, intent(inout) :: jgoto

private subroutine cnmn04(me, ii, xbar, eps, x1, y1, slope, x2, y2, x3, y3, x4, y4)

Routine to find first xbar>=eps corresponding to a minimum of a one-dimensional real function by polynomial interpolation.

Read more…

Arguments

Type IntentOptional Attributes Name
class(conmin_class), intent(inout) :: me
integer, intent(inout) :: ii

CALCULATION CONTROL:

Read more…
real(kind=wp), intent(out) :: xbar
real(kind=wp), intent(in) :: eps

may be negative

real(kind=wp), intent(in) :: x1
real(kind=wp), intent(in) :: y1
real(kind=wp), intent(in) :: slope
real(kind=wp), intent(in) :: x2
real(kind=wp), intent(in) :: y2
real(kind=wp), intent(in) :: x3
real(kind=wp), intent(in) :: y3
real(kind=wp), intent(in) :: x4
real(kind=wp), intent(in) :: y4

private subroutine cnmn05(me, g, df, a, s, b, c, slope, phi, isc, ic, ms1, nvc, n1, n2, n3, n4, n5)

Routine to solve direction finding problem in modified method of feasible directions.

Read more…

Arguments

Type IntentOptional Attributes Name
class(conmin_class), intent(inout) :: me
real(kind=wp), intent(inout) :: g(n2)
real(kind=wp), intent(inout) :: df(:)
real(kind=wp), intent(inout) :: a(n1,n3)
real(kind=wp), intent(inout) :: s(:)
real(kind=wp), intent(inout) :: b(n3,n3)
real(kind=wp), intent(inout) :: c(n4)
real(kind=wp), intent(inout) :: slope
real(kind=wp), intent(inout) :: phi
integer, intent(inout) :: isc(n2)
integer, intent(inout) :: ic(n3)
integer, intent(inout) :: ms1(n5)
integer, intent(inout) :: nvc
integer, intent(in) :: n1
integer, intent(in) :: n2
integer, intent(in) :: n3
integer, intent(in) :: n4
integer, intent(in) :: n5

private subroutine cnmn06(me, x, vlb, vub, g, scal, df, s, g1, g2, ctam, ctbm, slope, alp, a2, a3, a4, f1, f2, f3, cv1, cv2, cv3, cv4, alpca, alpfes, alpln, alpmin, alpnc, alpsav, alpsid, alptot, isc, ncal, nvc, icount, igood1, igood2, igood3, igood4, ibest, iii, nlnc, jgoto)

Routine to solve one-dimensional search problem for constrained function minimization.

Read more…

Arguments

Type IntentOptional Attributes Name
class(conmin_class), intent(inout) :: me
real(kind=wp), intent(inout) :: x(:)
real(kind=wp), intent(inout) :: vlb(:)
real(kind=wp), intent(inout) :: vub(:)
real(kind=wp), intent(inout) :: g(:)
real(kind=wp), intent(inout) :: scal(:)
real(kind=wp), intent(inout) :: df(:)
real(kind=wp), intent(inout) :: s(:)
real(kind=wp), intent(inout) :: g1(:)
real(kind=wp), intent(inout) :: g2(:)
real(kind=wp), intent(inout) :: ctam
real(kind=wp), intent(inout) :: ctbm
real(kind=wp), intent(inout) :: slope
real(kind=wp), intent(inout) :: alp
real(kind=wp), intent(inout) :: a2
real(kind=wp), intent(inout) :: a3
real(kind=wp), intent(inout) :: a4
real(kind=wp), intent(inout) :: f1
real(kind=wp), intent(inout) :: f2
real(kind=wp), intent(inout) :: f3
real(kind=wp), intent(inout) :: cv1
real(kind=wp), intent(inout) :: cv2
real(kind=wp), intent(inout) :: cv3
real(kind=wp), intent(inout) :: cv4
real(kind=wp), intent(inout) :: alpca
real(kind=wp), intent(inout) :: alpfes
real(kind=wp), intent(inout) :: alpln
real(kind=wp), intent(inout) :: alpmin
real(kind=wp), intent(inout) :: alpnc
real(kind=wp), intent(inout) :: alpsav
real(kind=wp), intent(inout) :: alpsid
real(kind=wp), intent(inout) :: alptot
integer, intent(inout) :: isc(:)
integer, intent(inout) :: ncal(2)
integer, intent(inout) :: nvc
integer, intent(inout) :: icount
integer, intent(inout) :: igood1
integer, intent(inout) :: igood2
integer, intent(inout) :: igood3
integer, intent(inout) :: igood4
integer, intent(inout) :: ibest
integer, intent(inout) :: iii
integer, intent(inout) :: nlnc
integer, intent(inout) :: jgoto

private subroutine cnmn07(me, ii, xbar, eps, x1, y1, x2, y2, x3, y3)

Routine to find first xbar>=epscorresponding to a real zero of a one-dimensional function by polynomial interpolation.

Read more…

Arguments

Type IntentOptional Attributes Name
class(conmin_class), intent(inout) :: me
integer, intent(inout) :: ii

CALCULATION CONTROL:

Read more…
real(kind=wp), intent(out) :: xbar
real(kind=wp), intent(in) :: eps

may be negative

real(kind=wp), intent(in) :: x1
real(kind=wp), intent(in) :: y1
real(kind=wp), intent(in) :: x2
real(kind=wp), intent(in) :: y2
real(kind=wp), intent(in) :: x3
real(kind=wp), intent(in) :: y3

private subroutine cnmn08(ndb, ner, c, ms1, b, n3, n4)

Routine to solve special linear problem for imposing s-transpose times s<=1 bounds in the modified method of feasible directions.

Read more…

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: ndb
integer, intent(out) :: ner

NER = ERROR FLAG. IF NER/=0 ON RETURN, PROCESS HAS NOT CONVERGED IN 5*NDB ITERATIONS.

real(kind=wp), intent(inout) :: c(n4)
integer, intent(out) :: ms1(:)

VECTOR MS1 IDENTIFIES THE SET OF BASIC VARIABLES.

real(kind=wp), intent(inout) :: b(n3,n3)
integer, intent(in) :: n3
integer, intent(in) :: n4