initialize the slsqp_solver class. see slsqp for more details.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(slsqp_solver), | intent(inout) | :: | me | |||
integer, | intent(in) | :: | n |
the number of variables, |
||
integer, | intent(in) | :: | m |
total number of constraints, |
||
integer, | intent(in) | :: | meq |
number of equality constraints, |
||
integer, | intent(in) | :: | max_iter |
maximum number of iterations |
||
real(kind=wp), | intent(in) | :: | acc |
accuracy |
||
procedure(func) | :: | f |
problem function |
|||
procedure(grad) | :: | g |
function to compute gradients (must be
associated if |
|||
real(kind=wp), | intent(in), | dimension(n) | :: | xl |
lower bounds on |
|
real(kind=wp), | intent(in), | dimension(n) | :: | xu |
upper bounds on |
|
logical, | intent(out) | :: | status_ok |
will be false if there were errors |
||
integer, | intent(in), | optional | :: | linesearch_mode |
1 = inexact (default), 2 = exact |
|
integer, | intent(in), | optional | :: | iprint |
unit number of status messages (default= |
|
procedure(iterfunc), | optional | :: | report |
user-defined procedure that will be called once per iteration |
||
real(kind=wp), | intent(in), | optional | :: | alphamin |
minimum alpha for linesearch [default 0.1] |
|
real(kind=wp), | intent(in), | optional | :: | alphamax |
maximum alpha for linesearch [default 1.0] |
|
integer, | intent(in), | optional | :: | gradient_mode |
how the gradients are to be computed:
Note that modes 1-3 do not respect the variable bounds. |
|
real(kind=wp), | intent(in), | optional | :: | gradient_delta |
perturbation step size (>epsilon) to compute the approximated
gradient by finite differences ( |
|
real(kind=wp), | intent(in), | optional | :: | tolf |
stopping criterion if then stop. |
|
real(kind=wp), | intent(in), | optional | :: | toldf |
stopping criterion if then stop |
|
real(kind=wp), | intent(in), | optional | :: | toldx |
stopping criterion if then stop |
|
integer, | intent(in), | optional | :: | max_iter_ls |
maximum number of iterations in the nnls problem |
|
integer, | intent(in), | optional | :: | nnls_mode |
Which NNLS method to use: |
|
real(kind=wp), | intent(in), | optional | :: | infinite_bound |
"infinity" for the upper and lower bounds.
if |
subroutine initialize_slsqp(me,n,m,meq,max_iter,acc,f,g,xl,xu,status_ok,& linesearch_mode,iprint,report,alphamin,alphamax,& gradient_mode,gradient_delta,tolf,toldf,toldx,& max_iter_ls,nnls_mode,infinite_bound) implicit none class(slsqp_solver),intent(inout) :: me integer,intent(in) :: n !! the number of variables, \( n \ge 1 \) integer,intent(in) :: m !! total number of constraints, \( m \ge 0 \) integer,intent(in) :: meq !! number of equality constraints, \( m_{eq} \ge 0 \) integer,intent(in) :: max_iter !! maximum number of iterations procedure(func) :: f !! problem function procedure(grad) :: g !! function to compute gradients (must be !! associated if `gradient_mode=0`) real(wp),dimension(n),intent(in) :: xl !! lower bounds on `x`. !! `xl(i)=NaN` (or `xl(i)<=-infinite_bound`) indicates to ignore `i`th bound real(wp),dimension(n),intent(in) :: xu !! upper bounds on `x`. !! `xu(i)=NaN` (or `xu(i)>=infinite_bound`) indicates to ignore `i`th bound real(wp),intent(in) :: acc !! accuracy logical,intent(out) :: status_ok !! will be false if there were errors integer,intent(in),optional :: linesearch_mode !! 1 = inexact (default), 2 = exact integer,intent(in),optional :: iprint !! unit number of status messages (default=`output_unit`) procedure(iterfunc),optional :: report !! user-defined procedure that will be called once per iteration real(wp),intent(in),optional :: alphamin !! minimum alpha for linesearch [default 0.1] real(wp),intent(in),optional :: alphamax !! maximum alpha for linesearch [default 1.0] integer,intent(in),optional :: gradient_mode !! how the gradients are to be computed: !! !! * 0 - use the user-supplied `g` subroutine. [default] !! * 1 - approximate by basic backward differences !! * 2 - approximate by basic forward differences !! * 3 - approximate by basic central differences !! !! Note that modes 1-3 do not respect the variable bounds. real(wp),intent(in),optional :: gradient_delta !! perturbation step size (>epsilon) to compute the approximated !! gradient by finite differences (`gradient_mode` 1-3). !! note that this is an absolute step that does not respect !! the `xl` or `xu` variable bounds. real(wp),intent(in),optional :: tolf !! stopping criterion if \( |f| < tolf \) then stop. real(wp),intent(in),optional :: toldf !! stopping criterion if \( |f_{n+1} - f_n| < toldf \) then stop real(wp),intent(in),optional :: toldx !! stopping criterion if \( ||x_{n+1} - x_n|| < toldx \) then stop integer,intent(in),optional :: max_iter_ls !! maximum number of iterations in the [[nnls]] problem integer,intent(in),optional :: nnls_mode !! Which NNLS method to use: !! !! 1. Use the original [[nnls]] !! 2. Use the newer [[bvls]] real(wp),intent(in),optional :: infinite_bound !! "infinity" for the upper and lower bounds. !! if `xl<=-infinite_bound` or `xu>=infinite_bound` !! then these bounds are considered nonexistant. !! If not present then `huge()` is used for this. integer :: n1,mineq,i status_ok = .false. call me%destroy() if (present(iprint)) me%iprint = iprint if (size(xl)/=size(xu) .or. size(xl)/=n) then call me%report_message('error: invalid upper or lower bound vector size') call me%report_message(' size(xl) =',ival=size(xl)) call me%report_message(' size(xu) =',ival=size(xu)) call me%report_message(' n =',ival=n) else if (meq<0 .or. meq>m) then call me%report_message('error: invalid meq value:', ival=meq) else if (m<0) then call me%report_message('error: invalid m value:', ival=m) else if (n<1) then call me%report_message('error: invalid n value:', ival=n) else if (any(xl>xu .and. .not. ieee_is_nan(xl) .and. .not. ieee_is_nan(xu))) then call me%report_message('error: lower bounds must be <= upper bounds.') do i=1,n if (xl(i)>xu(i) .and. .not. ieee_is_nan(xl(i)) .and. .not. ieee_is_nan(xu(i))) then call me%report_message(' xl(i)>xu(i) for variable',ival=i) end if end do else if (present(linesearch_mode)) then !two linesearch modes: select case (linesearch_mode) case(1) !inexact me%linesearch_mode = linesearch_mode case(2) !exact me%linesearch_mode = linesearch_mode case default call me%report_message('error: invalid linesearch_mode (must be 1 or 2): ',& ival=linesearch_mode) call me%destroy() return end select end if !optional linesearch bounds: if (present(alphamin)) me%alphamin = alphamin if (present(alphamax)) me%alphamax = alphamax !verify valid values for alphamin and alphamax: 0<alphamin<alphamax<=1 if (me%alphamin<=zero .or. me%alphamax<=zero .or. & me%alphamax<=me%alphamin .or. & me%alphamin>=one .or. me%alphamax>one) then call me%report_message('error: invalid values for alphamin or alphamax.') call me%report_message(' alphamin =',rval=me%alphamin) call me%report_message(' alphamax =',rval=me%alphamax) call me%destroy() return end if if (present(tolf)) me%tolf = tolf if (present(toldf)) me%toldf = toldf if (present(toldx)) me%toldx = toldx if (present(max_iter_ls)) me%max_iter_ls = max_iter_ls if (present(nnls_mode)) then select case (nnls_mode) case(1:2) me%nnls_mode = nnls_mode case default call me%report_message('error: invalid value for nnls_mode. defaulting to 1.') me%nnls_mode = 1 end select end if status_ok = .true. me%n = n me%m = m me%meq = meq me%max_iter = max_iter me%acc = acc me%f => f me%g => g if (present(report)) me%report => report allocate(me%xl(n)); me%xl = xl allocate(me%xu(n)); me%xu = xu !work arrays: n1 = n+1 mineq = m - meq + 2*n1 me%l_w = n1*(n1+1) + meq*(n1+1) + mineq*(n1+1) + & !for lsq (n1-meq+1)*(mineq+2) + 2*mineq + & !for lsi (n1+mineq)*(n1-meq) + 2*meq + n1 + & !for lsei n1*n/2 + 2*m + 3*n +3*n1 + 1 !for slsqpb allocate(me%w(me%l_w)) me%w = zero if (present(gradient_mode)) then me%gradient_mode = gradient_mode if (present(gradient_delta)) then me%gradient_delta = gradient_delta end if end if if (present(infinite_bound)) then me%infinite_bound = abs(infinite_bound) else me%infinite_bound = huge(one) end if end if end subroutine initialize_slsqp