This subroutine seeks the least value of a function of many variables, by a trust region method that forms quadratic models by interpolation.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| integer, | intent(in) | :: | n |
the number of variables and must be at least two |
||
| real(kind=wp), | intent(inout), | dimension(*) | :: | x |
Initial values of the variables must be set in X(1),X(2),...,X(N). They will be changed to the values that give the least calculated F. |
|
| real(kind=wp), | intent(in) | :: | rhobeg |
RHOBEG and RHOEND must be set to the initial and final values of a trust region radius, so both must be positive with RHOEND<=RHOBEG. Typically RHOBEG should be about one tenth of the greatest expected change to a variable, and RHOEND should indicate the accuracy that is required in the final values of the variables. |
||
| real(kind=wp), | intent(in) | :: | rhoend |
RHOBEG and RHOEND must be set to the initial and final values of a trust region radius, so both must be positive with RHOEND<=RHOBEG. Typically RHOBEG should be about one tenth of the greatest expected change to a variable, and RHOEND should indicate the accuracy that is required in the final values of the variables. |
||
| integer, | intent(in) | :: | iprint |
The value of IPRINT should be set to 0, 1, 2 or 3, which controls the amount of printing. Specifically, there is no output if IPRINT=0 and there is output only at the return if IPRINT=1. Otherwise, each new value of RHO is printed, with the best vector of variables so far and the corresponding value of the objective function. Further, each new value of F with its variables are output if IPRINT=3. |
||
| integer, | intent(in) | :: | maxfun |
upper bound on the number of calls of CALFUN. |
||
| procedure(func) | :: | calfun |
It must set F to the value of the objective function for the variables X(1),X(2),...,X(N). |
subroutine uobyqa (n, x, rhobeg, rhoend, iprint, maxfun, calfun) implicit none integer,intent(in) :: n !! the number of variables and must be at least two real(wp),intent(inout),dimension(*) :: x !! Initial values of the variables must be set in X(1),X(2),...,X(N). They !! will be changed to the values that give the least calculated F. real(wp),intent(in) :: rhobeg !! RHOBEG and RHOEND must be set to the initial and final values of a trust !! region radius, so both must be positive with RHOEND<=RHOBEG. Typically !! RHOBEG should be about one tenth of the greatest expected change to a !! variable, and RHOEND should indicate the accuracy that is required in !! the final values of the variables. real(wp),intent(in) :: rhoend !! RHOBEG and RHOEND must be set to the initial and final values of a trust !! region radius, so both must be positive with RHOEND<=RHOBEG. Typically !! RHOBEG should be about one tenth of the greatest expected change to a !! variable, and RHOEND should indicate the accuracy that is required in !! the final values of the variables. integer,intent(in) :: iprint !! The value of IPRINT should be set to 0, 1, 2 or 3, which controls the !! amount of printing. Specifically, there is no output if IPRINT=0 and !! there is output only at the return if IPRINT=1. Otherwise, each new !! value of RHO is printed, with the best vector of variables so far and !! the corresponding value of the objective function. Further, each new !! value of F with its variables are output if IPRINT=3. integer,intent(in) :: maxfun !! upper bound on the number of calls of CALFUN. procedure (func) :: calfun !! It must set F to the value of the objective !! function for the variables X(1),X(2),...,X(N). real(wp),dimension(:),allocatable :: w integer :: npt,ixb,ixo,ixn,ixp,ipq,ipl,ih,ig,id,ivl,iw ! Partition the working space array, so that different parts of it can be ! treated separately by the subroutine that performs the main calculation. npt = (n*n+3*n+2) / 2 ixb = 1 ixo = ixb + n ixn = ixo + n ixp = ixn + n ipq = ixp + n * npt ipl = ipq + npt - 1 ih = ipl + (npt-1) * npt ig = ih + n * n id = ig + n ivl = ih iw = id + n ! The array W will be used for working space allocate(w((N**4+8*N**3+23*N**2+42*N+max(2*N**2+4,18*N))/4)) call uobyqb (n, x, rhobeg, rhoend, iprint, maxfun, npt, w(ixb), w(ixo), w(ixn), & w(ixp), w(ipq), w(ipl), w(ih), w(ig), w(id), w(ivl), w(iw), calfun) deallocate(w) end subroutine uobyqa