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