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