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