CONOPT
Loading...
Searching...
No Matches
qp5.f90
Go to the documentation of this file.
1!> @file qp5.f90
2!! @ingroup FORT1THREAD_EXAMPLES
3!! @brief Similar to qp1 but uses 2nd derivatives computed using perturbations
4!!
5!! For more information about the individual callbacks, please have a look at the source code.
6
7! The current model is a simple QP model with a sparse Q matrix,
8! bounded variables, and one constraint. The number of superbasic
9! variables is larger than the default limit of 500 and the model
10! does not solve nicely. To solve it faster, there are these
11! possibilities:
12! 1. Increase the limit on the number of superbasics.
13! This can be done with a call to COIDEF_MaxSup or with an option.
14! Both approached are shown in extra solves in this file.
15! 2. Use directional 2nd derivatives, see qp2
16! 3. Use 2nd derivatives as a matrix, see qp3
17! 4. 2 and 3 combined, see qp4
18! 5. Use 2nd derivatives computed using perturbations, see qp5
19!
20
21!> Main program. A simple setup and call of CONOPT for a QP model
22!!
23Program qp
24
25 Use proginfo
26 Use coidef
27 Use qp_data
28 implicit None
29!
30! Declare the user callback routines as Integer, External:
31!
32 Integer, External :: qp_readmatrix ! Mandatory Matrix definition routine defined below
33 Integer, External :: qp_fdeval ! Function and Derivative evaluation routine
34 ! needed a nonlinear model.
35 Integer, External :: std_status ! Standard callback for displaying solution status
36 Integer, External :: std_solution ! Standard callback for displaying solution values
37 Integer, External :: std_message ! Standard callback for managing messages
38 Integer, External :: std_errmsg ! Standard callback for managing error messages
39#if defined(itl)
40!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: QP_ReadMatrix
41!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: QP_FDEval
42!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Status
43!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Solution
44!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Message
45!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_ErrMsg
46#endif
47!
48! Control vector
49!
50 INTEGER :: numcallback
51 INTEGER, Dimension(:), Pointer :: cntvect
52 INTEGER :: coi_error
53!
54! Local variables used to define Q
55!
56 Integer :: i,j
57
58 Call startup
59!
60! Create and initialize a Control Vector
61!
62 numcallback = coidef_size()
63 Allocate( cntvect(numcallback) )
64 coi_error = coidef_inifort( cntvect )
65!
66! Tell CONOPT about the size of the model by populating the Control Vector:
67! We will create a QP model with one constraint, sum(i, x(i) ) = 1
68! and NN variables. NN and the Q matrix are declared in Module QP_Data.
69!
70 j = 1
71 do i = 1, nn
72 q(j) = 1.d0
73 qrow(j) = i
74 qcol(j) = i
75 j = j + 1
76 if ( i < nn ) then
77 q(j) = 0.1d0
78 qrow(j) = i+1
79 qcol(j) = i
80 j = j + 1
81 endif
82 enddo
83 do i = 1, nn
84 target(i) = 10.d0
85 enddo
86
87 coi_error = max( coi_error, coidef_numvar( cntvect, nn ) ) ! NN variables
88 coi_error = max( coi_error, coidef_numcon( cntvect, 2 ) ) ! 1 constraint + 1 objective
89 coi_error = max( coi_error, coidef_numnz( cntvect, 2*nn ) ) ! 2*NN nonzeros in the Jacobian
90 coi_error = max( coi_error, coidef_numnlnz( cntvect, nn ) ) ! NN of which are nonlinear
91 coi_error = max( coi_error, coidef_optdir( cntvect, -1 ) ) ! Minimize
92 coi_error = max( coi_error, coidef_objcon( cntvect, 1 ) ) ! Objective is constraint 1
93 coi_error = max( coi_error, coidef_optfile( cntvect, 'qp5.opt' ) )
94!
95! Tell CONOPT about the callback routines:
96!
97 coi_error = max( coi_error, coidef_readmatrix( cntvect, qp_readmatrix ) )
98 coi_error = max( coi_error, coidef_fdeval( cntvect, qp_fdeval ) )
99 coi_error = max( coi_error, coidef_status( cntvect, std_status ) )
100 coi_error = max( coi_error, coidef_solution( cntvect, std_solution ) )
101 coi_error = max( coi_error, coidef_message( cntvect, std_message ) )
102 coi_error = max( coi_error, coidef_errmsg( cntvect, std_errmsg ) )
103
104#if defined(LICENSE_INT_1) && defined(LICENSE_INT_2) && defined(LICENSE_INT_3) && defined(LICENSE_TEXT)
105 coi_error = max( coi_error, coidef_license( cntvect, license_int_1, license_int_2, license_int_3, license_text) )
106#endif
107
108 If ( coi_error .ne. 0 ) THEN
109 write(*,*)
110 write(*,*) '**** Fatal Error while loading CONOPT Callback routines.'
111 write(*,*)
112 call flog( "Skipping Solve due to setup errors", 1 )
113 ENDIF
114!
115! Start CONOPT:
116!
117 coi_error = coi_solve( cntvect )
118
119 write(*,*)
120 write(*,*) 'End of QP example 5. Return code=',coi_error
121
122 If ( coi_error /= 0 ) then
123 call flog( "Errors encountered during solution", 1 )
124 elseif ( stacalls == 0 .or. solcalls == 0 ) then
125 call flog( "Status or Solution routine was not called", 1 )
126 elseif ( sstat /= 1 .or. mstat /= 2 ) then
127 call flog( "Solver and Model Status was not as expected (1,2)", 1 )
128 endif
129
130 if ( coi_free(cntvect) /= 0 ) call flog( "Error while freeing control vector",1)
131
132 call flog( "Successful Solve", 0 )
133
134End Program qp
135
136!> Define information about the model
137!!
138!! @include{doc} readMatrix_params.dox
139Integer Function qp_readmatrix( lower, curr, upper, vsta, type, rhs, esta, &
140 colsta, rowno, value, nlflag, n, m, nz, &
141 usrmem )
142#if defined(itl)
143!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: QP_ReadMatrix
144#endif
145 Use qp_data
146 implicit none
147 integer, intent (in) :: n ! number of variables
148 integer, intent (in) :: m ! number of constraints
149 integer, intent (in) :: nz ! number of nonzeros
150 real*8, intent (in out), dimension(n) :: lower ! vector of lower bounds
151 real*8, intent (in out), dimension(n) :: curr ! vector of initial values
152 real*8, intent (in out), dimension(n) :: upper ! vector of upper bounds
153 integer, intent (in out), dimension(n) :: vsta ! vector of initial variable status
154 ! (not defined here)
155 integer, intent (out), dimension(m) :: type ! vector of equation types
156 integer, intent (in out), dimension(m) :: esta ! vector of initial equation status
157 ! (not defined here)
158 real*8, intent (in out), dimension(m) :: rhs ! vector of right hand sides
159 integer, intent (in out), dimension(n+1) :: colsta ! vector with start of column indices
160 integer, intent (out), dimension(nz) :: rowno ! vector of row numbers
161 integer, intent (in out), dimension(nz) :: nlflag ! vector of nonlinearity flags
162 real*8, intent (in out), dimension(nz) :: value ! vector of matrix values
163 real*8 usrmem(*) ! optional user memory
164 integer :: i, j
165!
166! Information about Variables:
167! Default: Lower = -Inf, Curr = 0, and Upper = +inf.
168! Default: the status information in Vsta is not used.
169!
170! Nondefault: Lower bound = 0 for all variables
171!
172 do i = 1, nn
173 lower(i) = 0.0d0
174 enddo
175!
176! Information about Constraints:
177! Default: Rhs = 0
178! Default: the status information in Esta and the function
179! value in FV are not used.
180! Default: Type: There is no default.
181! 0 = Equality,
182! 1 = Greater than or equal,
183! 2 = Less than or equal,
184! 3 = Non binding.
185!
186! Constraint 1 (Objective)
187! Rhs = 0.0 and type Non binding
188!
189 type(1) = 3
190!
191! Constraint 2 (Sum of all vars = 1)
192! Rhs = 1 and type Equality
193!
194 rhs(2) = 1.d0
195 type(2) = 0
196!
197! Information about the Jacobian. We use the standard method with
198! Rowno, Value, Nlflag and Colsta and we do not use Colno.
199!
200! Colsta = Start of column indices (No Defaults):
201! Rowno = Row indices
202! Value = Value of derivative (by default only linear
203! derivatives are used)
204! Nlflag = 0 for linear and 1 for nonlinear derivative
205! (not needed for completely linear models)
206!
207 j = 1 ! counter for current nonzero
208 do i = 1, nn
209 colsta(i) = j
210 rowno(j) = 1
211 nlflag(j) = 1
212 j = j + 1
213 rowno(j) = 2
214 value(j) = 1.d0
215 nlflag(j) = 0
216 j = j + 1
217 enddo
218 colsta(nn+1) = j
219
220 qp_readmatrix = 0 ! Return value means OK
221
222end Function qp_readmatrix
223
224!> Compute nonlinear terms and non-constant Jacobian elements
225!!
226!! @include{doc} fdeval_params.dox
227Integer Function qp_fdeval( x, g, jac, rowno, jcnm, mode, ignerr, errcnt, &
228 n, nz, thread, usrmem )
229#if defined(itl)
230!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: QP_FDEval
231#endif
232 Use qp_data
233 implicit none
234 integer, intent (in) :: n ! number of variables
235 integer, intent (in) :: rowno ! number of the row to be evaluated
236 integer, intent (in) :: nz ! number of nonzeros in this row
237 real*8, intent (in), dimension(n) :: x ! vector of current solution values
238 real*8, intent (in out) :: g ! constraint value
239 real*8, intent (in out), dimension(n) :: jac ! vector of derivatives for current constraint
240 integer, intent (in), dimension(nz) :: jcnm ! list of variables that appear nonlinearly
241 ! in this row. Ffor information only.
242 integer, intent (in) :: mode ! evaluation mode: 1 = function value
243 ! 2 = derivatives, 3 = both
244 integer, intent (in) :: ignerr ! if 1 then errors can be ignored as long
245 ! as errcnt is incremented
246 integer, intent (in out) :: errcnt ! error counter to be incremented in case
247 ! of function evaluation errors.
248 integer, intent (in) :: thread
249 real*8 usrmem(*) ! optional user memory
250 Integer :: i, j, k
251 real*8 :: sum
252!
253! Row 1: The objective function is nonlinear
254!
255 if ( rowno .eq. 1 ) then
256!
257! Mode = 1 or 3. Function value: G = x * Q * x / 2
258!
259 if ( mode .eq. 1 .or. mode .eq. 3 ) then
260 sum = 0.d0
261 do k = 1, nq
262 i = qrow(k)
263 j = qcol(k)
264 if ( i == j ) then
265 sum = sum + (x(i)-target(i))*q(k)*(x(j)-target(j))
266 else
267 sum = sum + 2*(x(i)-target(i))*q(k)*(x(j)-target(j))
268 endif
269 enddo
270 g = sum / 2.d0
271 endif
272!
273! Mode = 2 or 3: Derivative values: Q*x
274!
275 if ( mode .eq. 2 .or. mode .eq. 3 ) then
276 do i = 1, nn
277 jac(i) = 0
278 enddo
279 do k = 1, nq
280 i = qrow(k)
281 j = qcol(k)
282 if ( i == j ) then
283 jac(i) = jac(i) + q(k) * (x(i)-target(i))
284 else
285 jac(i) = jac(i) + q(k) * (x(j)-target(j))
286 jac(j) = jac(j) + q(k) * (x(i)-target(i))
287 endif
288 enddo
289 endif
290!
291! Row = 2: The row is linear and will not be called.
292!
293 endif
294 qp_fdeval = 0
295
296end Function qp_fdeval
integer function std_solution(xval, xmar, xbas, xsta, yval, ymar, ybas, ysta, n, m, usrmem)
Definition comdecl.f90:128
integer function std_status(modsta, solsta, iter, objval, usrmem)
Definition comdecl.f90:82
integer function std_message(smsg, dmsg, nmsg, llen, usrmem, msgv)
Definition comdecl.f90:203
integer function std_errmsg(rowno, colno, posno, msglen, usrmem, msg)
Definition comdecl.f90:248
integer function coidef_fdeval(cntvect, coi_fdeval)
define callback routine for performing function and derivative evaluations.
integer function coidef_errmsg(cntvect, coi_errmsg)
define callback routine for returning error messages for row, column or Jacobian elements.
integer function coidef_message(cntvect, coi_message)
define callback routine for handling messages returned during the solution process.
integer function coidef_readmatrix(cntvect, coi_readmatrix)
define callback routine for providing the matrix data to CONOPT.
integer function coidef_status(cntvect, coi_status)
define callback routine for returning the completion status.
integer function coidef_solution(cntvect, coi_solution)
define callback routine for returning the final solution values.
integer function coidef_optfile(cntvect, optfile)
define callback routine for defining an options file.
integer function coidef_license(cntvect, licint1, licint2, licint3, licstring)
define the License Information.
Definition coistart.f90:680
integer function coidef_numvar(cntvect, numvar)
defines the number of variables in the model.
Definition coistart.f90:358
integer function coidef_objcon(cntvect, objcon)
defines the Objective Constraint.
Definition coistart.f90:629
integer function coidef_numnz(cntvect, numnz)
defines the number of nonzero elements in the Jacobian.
Definition coistart.f90:437
integer function coidef_optdir(cntvect, optdir)
defines the Optimization Direction.
Definition coistart.f90:552
integer function coidef_numnlnz(cntvect, numnlnz)
defines the Number of Nonlinear Nonzeros.
Definition coistart.f90:476
integer function coidef_numcon(cntvect, numcon)
defines the number of constraints in the model.
Definition coistart.f90:398
integer function coidef_size()
returns the size the Control Vector must have, measured in standard Integer units.
Definition coistart.f90:176
integer function coidef_inifort(cntvect)
initialisation method for Fortran applications.
Definition coistart.f90:314
integer function coi_solve(cntvect)
method for starting the solving process of CONOPT.
Definition coistart.f90:14
integer solcalls
Definition comdecl.f90:9
integer sstat
Definition comdecl.f90:12
integer stacalls
Definition comdecl.f90:8
subroutine flog(msg, code)
Definition comdecl.f90:56
integer mstat
Definition comdecl.f90:11
subroutine startup
Definition comdecl.f90:35
real(8) target
integer, parameter nn
real(8), dimension(nn, nn) q
program qp
Main program. A simple setup and call of CONOPT for a QP model.
Definition qp1.f90:26
integer function qp_readmatrix(lower, curr, upper, vsta, type, rhs, esta, colsta, rowno, value, nlflag, n, m, nz, usrmem)
Define information about the model.
Definition qp1.f90:209
integer function qp_fdeval(x, g, jac, rowno, jcnm, mode, ignerr, errcnt, n, nz, thread, usrmem)
Compute nonlinear terms and non-constant Jacobian elements.
Definition qp1.f90:296