CONOPT
Loading...
Searching...
No Matches
qp3.f90
Go to the documentation of this file.
1!> @file qp3.f90
2!! @ingroup FORT1THREAD_EXAMPLES
3!! @brief Similar to qp1 but uses 2nd derivatives as a matrix
4!!
5!! For more information about the individual callbacks, please have a look at the source code.
6
7#if defined(_WIN32) && !defined(_WIN64)
8#define dec_directives_win32
9#endif
10
11!> Main program. A simple setup and call of CONOPT for a QP model
12!!
13Program qp
14
16 Use conopt
17 Use qp_data
18 implicit None
19!
20! Declare the user callback routines as Integer, External:
21!
22 Integer, External :: qp_readmatrix ! Mandatory Matrix definition routine defined below
23 Integer, External :: qp_fdeval ! Function and Derivative evaluation routine
24 ! needed a nonlinear model.
25 Integer, External :: qp_2dlagrstr ! 2nd derivative routine for structure
26 Integer, External :: qp_2dlagrval ! 2nd derivative routine for values
27 Integer, External :: std_status ! Standard callback for displaying solution status
28 Integer, External :: std_solution ! Standard callback for displaying solution values
29 Integer, External :: std_message ! Standard callback for managing messages
30 Integer, External :: std_errmsg ! Standard callback for managing error messages
31#ifdef dec_directives_win32
32!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: QP_ReadMatrix
33!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: QP_FDEval
34!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: QP_2DLagrStr
35!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: QP_2DLagrVal
36!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Status
37!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Solution
38!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Message
39!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_ErrMsg
40#endif
41!
42! Control vector
43!
44 INTEGER, Dimension(:), Pointer :: cntvect
45 INTEGER :: coi_error
46!
47! Local variables used to define Q
48!
49 Integer :: i,j
50
51 call startup
52!
53! Create and initialize a Control Vector
54!
55 coi_error = coi_create( cntvect )
56!
57! Tell CONOPT about the size of the model by populating the Control Vector:
58! We will create a QP model with one constraint, sum(i, x(i) ) = 1
59! and NN variables. NN and the Q matrix are defined in Module QP_Data.
60!
61 j = 1
62 do i = 1, nn
63 q(j) = 1.d0
64 qrow(j) = i
65 qcol(j) = i
66 j = j + 1
67 if ( i < nn ) then
68 q(j) = 0.1d0
69 qrow(j) = i+1
70 qcol(j) = i
71 j = j + 1
72 endif
73 enddo
74 do i = 1, nn
75 target(i) = 10.d0
76 enddo
77
78 coi_error = max( coi_error, coidef_numvar( cntvect, nn ) ) ! NN variables
79 coi_error = max( coi_error, coidef_numcon( cntvect, 2 ) ) ! 1 constraint + 1 objective
80 coi_error = max( coi_error, coidef_numnz( cntvect, 2*nn ) ) ! 2*NN nonzeros in the Jacobian
81 coi_error = max( coi_error, coidef_numnlnz( cntvect, nn ) ) ! NN of which are nonlinear
82 coi_error = max( coi_error, coidef_numhess( cntvect, nq ) ) ! NQ elements in the Hessian
83 coi_error = max( coi_error, coidef_optdir( cntvect, -1 ) ) ! Minimize
84 coi_error = max( coi_error, coidef_objcon( cntvect, 1 ) ) ! Objective is constraint 1
85!
86! Tell CONOPT about the callback routines:
87!
88 coi_error = max( coi_error, coidef_readmatrix( cntvect, qp_readmatrix ) )
89 coi_error = max( coi_error, coidef_fdeval( cntvect, qp_fdeval ) )
90 coi_error = max( coi_error, coidef_2dlagrstr( cntvect, qp_2dlagrstr ) )
91 coi_error = max( coi_error, coidef_2dlagrval( cntvect, qp_2dlagrval ) )
92 coi_error = max( coi_error, coidef_status( cntvect, std_status ) )
93 coi_error = max( coi_error, coidef_solution( cntvect, std_solution ) )
94 coi_error = max( coi_error, coidef_message( cntvect, std_message ) )
95 coi_error = max( coi_error, coidef_errmsg( cntvect, std_errmsg ) )
96
97#if defined(CONOPT_LICENSE_INT_1) && defined(CONOPT_LICENSE_INT_2) && defined(CONOPT_LICENSE_INT_3) && defined(CONOPT_LICENSE_TEXT)
98 coi_error = max( coi_error, coidef_license( cntvect, conopt_license_int_1, conopt_license_int_2, conopt_license_int_3, conopt_license_text) )
99#endif
100
101 If ( coi_error .ne. 0 ) THEN
102 write(*,*)
103 write(*,*) '**** Fatal Error while loading CONOPT Callback routines.'
104 write(*,*)
105 call flog( "Skipping Solve due to setup errors", 1 )
106 ENDIF
107!
108! Start CONOPT:
109!
110 coi_error = coi_solve( cntvect )
111
112 write(*,*)
113 write(*,*) 'End of QP example 3. Return code=',coi_error
114
115 If ( coi_error /= 0 ) then
116 call flog( "Errors encountered during solution", 1 )
117 elseif ( stacalls == 0 .or. solcalls == 0 ) then
118 call flog( "Status or Solution routine was not called", 1 )
119 elseif ( sstat /= 1 .or. mstat /= 2 ) then
120 call flog( "Solver and Model Status was not as expected (1,2)", 1 )
121 endif
122
123 if ( coi_free(cntvect) /= 0 ) call flog( "Error while freeing control vector",1)
124
125 call flog( "Successful Solve", 0 )
126!
127! Free solution memory
128!
129 call finalize
130
131End Program qp
132
133!> Define information about the model
134!!
135!! @include{doc} readMatrix_params.dox
136Integer Function qp_readmatrix( lower, curr, upper, vsta, type, rhs, esta, &
137 colsta, rowno, value, nlflag, n, m, nz, &
138 usrmem )
139#ifdef dec_directives_win32
140!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: QP_ReadMatrix
141#endif
142 Use qp_data
143 implicit none
144 integer, intent (in) :: n ! number of variables
145 integer, intent (in) :: m ! number of constraints
146 integer, intent (in) :: nz ! number of nonzeros
147 real*8, intent (in out), dimension(n) :: lower ! vector of lower bounds
148 real*8, intent (in out), dimension(n) :: curr ! vector of initial values
149 real*8, intent (in out), dimension(n) :: upper ! vector of upper bounds
150 integer, intent (in out), dimension(n) :: vsta ! vector of initial variable status
151 ! (not defined here)
152 integer, intent (out), dimension(m) :: type ! vector of equation types
153 integer, intent (in out), dimension(m) :: esta ! vector of initial equation status
154 ! (not defined here)
155 real*8, intent (in out), dimension(m) :: rhs ! vector of right hand sides
156 integer, intent (in out), dimension(n+1) :: colsta ! vector with start of column indices
157 integer, intent (out), dimension(nz) :: rowno ! vector of row numbers
158 integer, intent (in out), dimension(nz) :: nlflag ! vector of nonlinearity flags
159 real*8, intent (in out), dimension(nz) :: value ! vector of matrix values
160 real*8 usrmem(*) ! optional user memory
161 integer :: i, j
162!
163! Information about Variables:
164! Default: Lower = -Inf, Curr = 0, and Upper = +inf.
165! Default: the status information in Vsta is not used.
166!
167! Nondefault: Lower bound = 0 for all variables
168!
169 do i = 1, nn
170 lower(i) = 0.0d0
171 enddo
172!
173! Information about Constraints:
174! Default: Rhs = 0
175! Default: the status information in Esta and the function
176! value in FV are not used.
177! Default: Type: There is no default.
178! 0 = Equality,
179! 1 = Greater than or equal,
180! 2 = Less than or equal,
181! 3 = Non binding.
182!
183! Constraint 1 (Objective)
184! Rhs = 0.0 and type Non binding
185!
186 type(1) = 3
187!
188! Constraint 2 (Sum of all vars = 1)
189! Rhs = 1 and type Equality
190!
191 rhs(2) = 1.d0
192 type(2) = 0
193!
194! Information about the Jacobian. CONOPT expects a columnwise
195! representation in Rowno, Value, Nlflag and Colsta.
196!
197! Colsta = Start of column indices (No Defaults):
198! Rowno = Row indices
199! Value = Value of derivative (by default only linear
200! derivatives are used)
201! Nlflag = 0 for linear and 1 for nonlinear derivative
202! (not needed for completely linear models)
203!
204 j = 1 ! counter for current nonzero
205 do i = 1, nn
206 colsta(i) = j
207 rowno(j) = 1
208 nlflag(j) = 1
209 j = j + 1
210 rowno(j) = 2
211 value(j) = 1.d0
212 nlflag(j) = 0
213 j = j + 1
214 enddo
215 colsta(nn+1) = j
216
217 qp_readmatrix = 0 ! Return value means OK
218
219end Function qp_readmatrix
220
221!> Compute nonlinear terms and non-constant Jacobian elements
222!!
223!! @include{doc} fdeval_params.dox
224Integer Function qp_fdeval( x, g, jac, rowno, jcnm, mode, ignerr, errcnt, &
225 n, nz, thread, usrmem )
226#ifdef dec_directives_win32
227!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: QP_FDEval
228#endif
229 Use qp_data
230 implicit none
231 integer, intent (in) :: n ! number of variables
232 integer, intent (in) :: rowno ! number of the row to be evaluated
233 integer, intent (in) :: nz ! number of nonzeros in this row
234 real*8, intent (in), dimension(n) :: x ! vector of current solution values
235 real*8, intent (in out) :: g ! constraint value
236 real*8, intent (in out), dimension(n) :: jac ! vector of derivatives for current constraint
237 integer, intent (in), dimension(nz) :: jcnm ! list of variables that appear nonlinearly
238 ! in this row. Ffor information only.
239 integer, intent (in) :: mode ! evaluation mode: 1 = function value
240 ! 2 = derivatives, 3 = both
241 integer, intent (in) :: ignerr ! if 1 then errors can be ignored as long
242 ! as errcnt is incremented
243 integer, intent (in out) :: errcnt ! error counter to be incremented in case
244 ! of function evaluation errors.
245 integer, intent (in) :: thread
246 real*8 usrmem(*) ! optional user memory
247 Integer :: i,j, k
248 real*8 :: sum
249!
250! Row 1: The objective function is nonlinear
251!
252 if ( rowno .eq. 1 ) then
253!
254! Mode = 1 or 3. Function value: G = x * Q * x / 2
255!
256 if ( mode .eq. 1 .or. mode .eq. 3 ) then
257 sum = 0.d0
258 do k = 1, nq
259 i = qrow(k)
260 j = qcol(k)
261 if ( i == j ) then
262 sum = sum + (x(i)-target(i))*q(k)*(x(j)-target(j))
263 else
264 sum = sum + 2*(x(i)-target(i))*q(k)*(x(j)-target(j))
265 endif
266 enddo
267 g = sum / 2.d0
268 endif
269!
270! Mode = 2 or 3: Derivative values: Q*x
271!
272 if ( mode .eq. 2 .or. mode .eq. 3 ) then
273 do i = 1, nn
274 jac(i) = 0
275 enddo
276 do k = 1, nq
277 i = qrow(k)
278 j = qcol(k)
279 if ( i == j ) then
280 jac(i) = jac(i) + q(k) * (x(i)-target(i))
281 else
282 jac(i) = jac(i) + q(k) * (x(j)-target(j))
283 jac(j) = jac(j) + q(k) * (x(i)-target(i))
284 endif
285 enddo
286 endif
287!
288! Row = 2: The row is linear and will not be called.
289!
290 endif
291 qp_fdeval = 0
292
293end Function qp_fdeval
294
295!> Specify the structure of the Lagrangian of the Hessian
296!!
297!! @include{doc} 2DLagrStr_params.dox
298Integer Function qp_2dlagrstr( ROWNO, COLNO, NODRV, N, M, NHESS, USRMEM )
299#ifdef dec_directives_win32
300!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: QP_2DLagrStr
301#endif
302 Use qp_data
303 implicit none
304 INTEGER n, m, nhess, nodrv
305 real*8 usrmem(*)
306 INTEGER rowno(nhess), colno(nhess)
307
308 Integer :: j
309
310 Do j = 1, nq
311 rowno(j) = qrow(j)
312 colno(j) = qcol(j)
313 Enddo
314 qp_2dlagrstr = 0
315
316End function qp_2dlagrstr
317
318!> Compute the Lagrangian of the Hessian
319!!
320!! @include{doc} 2DLagrVal_params.dox
321Integer Function qp_2dlagrval( X, U, ROWNO, COLNO, VALUE, NODRV, N, M, NHESS, USRMEM )
322#ifdef dec_directives_win32
323!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: QP_2DLagrVal
324#endif
325 Use qp_data
326 implicit none
327 INTEGER n, m, nhess, nodrv
328 real*8 x(n), u(m), value(nhess), usrmem(*)
329 INTEGER rowno(nhess), colno(nhess)
330
331 Integer :: j
332
333 Do j = 1, nq
334 value(j) = q(j)*u(1)
335 Enddo
336 qp_2dlagrval = 0
337
338End function qp_2dlagrval
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_2dlagrstr(cntvect, coi_2dlagrstr)
define callback routine for providing the structure of the second derivatives of the Lagrangian.
Definition conopt.f90:1547
integer(c_int) function coidef_2dlagrval(cntvect, coi_2dlagrval)
define callback routine for computing the values of the second derivatives of the Lagrangian.
Definition conopt.f90:1573
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_numhess(cntvect, numhess)
defines the Number of Hessian Nonzeros.
Definition conopt.f90:188
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
integer function qp_2dlagrval(x, u, rowno, colno, value, nodrv, n, m, nhess, usrmem)
Compute the Lagrangian of the Hessian.
Definition qp3.f90:302
integer function qp_2dlagrstr(rowno, colno, nodrv, n, m, nhess, usrmem)
Specify the structure of the Lagrangian of the Hessian.
Definition qp3.f90:282