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