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