CONOPT
Loading...
Searching...
No Matches
sqp1.f90
Go to the documentation of this file.
1!> @file sqp1.f90
2!! @ingroup FORT1THREAD_EXAMPLES
3!!
4!!
5!! SQP1: 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!> Main program. A simple setup and call of CONOPT for a QP model
12!!
13Program qp
14
15 Use proginfo
16 Use coidef
17 Use qp_data_sp
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#if defined(itl)
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 :: numcallback
45 INTEGER, Dimension(:), Pointer :: cntvect
46 INTEGER :: coi_error
47!
48! Local variables used to define Q
49!
50 Integer :: i,j,k
51
52 call startup
53!
54! We will create a QP model with one constraint, sum(i, x(i) ) = 1
55! and NN variables. NN and the Q matrix are defined in Module QP_Data_Sp.
56!
57 j = 0
58 do i = 1, nn
59 j = j + 1
60 q(j) = 1.d0
61 qrow(j) = i
62 qcol(j) = i
63 do k = i+1, nn, 2
64 j = j + 1
65 q(j) = 0.5d0/nn ! Make sure the matrix is diagonally dominant.
66 qrow(j) = k
67 qcol(j) = i
68 enddo
69 enddo
70 nq = j
71 do i = 1, nn
72 target(i) = 0.1d0
73 enddo
74!
75! Create and initialize a Control Vector
76!
77 numcallback = coidef_size()
78 Allocate( cntvect(numcallback) )
79 coi_error = coidef_inifort( cntvect )
80!
81! Tell CONOPT about the size of the model by populating the Control Vector:
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, 1000.d0 ) ) ! Allow a dense Hessian in 2DLagr
91 coi_error = max( coi_error, coidef_optfile( cntvect, 'sqp1.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 write(*,*)
120 write(*,*) 'End of SQP example 1 case 1. Return code=',coi_error
121
122 If ( coi_error /= 0 ) then
123 call flog( "Errors encountered during solution of case 1", 1 )
124 elseif ( stacalls == 0 .or. solcalls == 0 ) then
125 call flog( "Status or Solution routine was not called in case 1", 1 )
126 elseif ( sstat /= 1 .or. mstat /= 2 ) then
127 call flog( "Solver and Model Status was not as expected (1,2) in case 1", 1 )
128 endif
129!
130! In case 2 we increate the targets. The objective and derivatives will become larger.
131
132 do i = 1, nn
133 target(i) = 1.0d0
134 enddo
135 stacalls = 0; solcalls = 0; sstat = 0; mstat = 0
136!
137! Call CONOPT again
138!
139 coi_error = coi_solve( cntvect )
140
141 write(*,*)
142 write(*,*) 'End of SQP example 1 case 2. Return code=',coi_error
143
144 If ( coi_error /= 0 ) then
145 call flog( "Errors encountered during solution of case 2", 1 )
146 elseif ( stacalls == 0 .or. solcalls == 0 ) then
147 call flog( "Status or Solution routine was not called in case 2", 1 )
148 elseif ( sstat /= 1 .or. mstat /= 2 ) then
149 call flog( "Solver and Model Status was not as expected (1,2) in case 2", 1 )
150 endif
151!
152! In case 3 we increate the targets again. The objective and derivatives will become much larger.
153
154 do i = 1, nn
155 target(i) = 50.0d0
156 enddo
157 stacalls = 0; solcalls = 0; sstat = 0; mstat = 0
158!
159! Call CONOPT again
160!
161 coi_error = coi_solve( cntvect )
162
163 write(*,*)
164 write(*,*) 'End of SQP example 1 case 3. Return code=',coi_error
165
166 If ( coi_error /= 0 ) then
167 call flog( "Errors encountered during solution of case 3", 1 )
168 elseif ( stacalls == 0 .or. solcalls == 0 ) then
169 call flog( "Status or Solution routine was not called in case 3", 1 )
170 elseif ( sstat /= 1 .or. mstat /= 2 ) then
171 call flog( "Solver and Model Status was not as expected (1,2) in case 3", 1 )
172 endif
173
174 if ( coi_free(cntvect) /= 0 ) call flog( "Error while freeing control vector",1)
175
176 call flog( "Successful Solve", 0 )
177
178End Program qp
179!
180! ============================================================================
181! Define information about the model:
182!
183
184!> Define information about the model
185!!
186!! @include{doc} readMatrix_params.dox
187Integer Function qp_readmatrix( lower, curr, upper, vsta, type, rhs, esta, &
188 colsta, rowno, value, nlflag, n, m, nz, &
189 usrmem )
190#if defined(itl)
191!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: QP_ReadMatrix
192#endif
193 Use qp_data_sp
194 implicit none
195 integer, intent (in) :: n ! number of variables
196 integer, intent (in) :: m ! number of constraints
197 integer, intent (in) :: nz ! number of nonzeros
198 real*8, intent (in out), dimension(n) :: lower ! vector of lower bounds
199 real*8, intent (in out), dimension(n) :: curr ! vector of initial values
200 real*8, intent (in out), dimension(n) :: upper ! vector of upper bounds
201 integer, intent (in out), dimension(n) :: vsta ! vector of initial variable status
202 ! (not defined here)
203 integer, intent (out), dimension(m) :: type ! vector of equation types
204 integer, intent (in out), dimension(m) :: esta ! vector of initial equation status
205 ! (not defined here)
206 real*8, intent (in out), dimension(m) :: rhs ! vector of right hand sides
207 integer, intent (in out), dimension(n+1) :: colsta ! vector with start of column indices
208 integer, intent (out), dimension(nz) :: rowno ! vector of row numbers
209 integer, intent (in out), dimension(nz) :: nlflag ! vector of nonlinearity flags
210 real*8, intent (in out), dimension(nz) :: value ! vector of matrix values
211 real*8 usrmem(*) ! optional user memory
212 integer :: i, j
213!
214! Information about Variables:
215! Default: Lower = -Inf, Curr = 0, and Upper = +inf.
216! Default: the status information in Vsta is not used.
217!
218! Nondefault: Lower bound = 0 for all variables
219!
220 do i = 1, nn
221 lower(i) = 0.0d0
222 enddo
223!
224! Information about Constraints:
225! Default: Rhs = 0
226! Default: the status information in Esta and the function
227! value in FV are not used.
228! Default: Type: There is no default.
229! 0 = Equality,
230! 1 = Greater than or equal,
231! 2 = Less than or equal,
232! 3 = Non binding.
233!
234! Constraint 1 (Objective)
235! Rhs = 0.0 and type Non binding
236!
237 type(1) = 3
238!
239! Constraint 2 (Sum of all vars = 1)
240! Rhs = 1 and type Equality
241!
242 rhs(2) = 1.d0
243 type(2) = 0
244!
245! Information about the Jacobian. We use the standard method with
246! Rowno, Value, Nlflag and Colsta and we do not use Colno.
247!
248! Colsta = Start of column indices (No Defaults):
249! Rowno = Row indices
250! Value = Value of derivative (by default only linear
251! derivatives are used)
252! Nlflag = 0 for linear and 1 for nonlinear derivative
253! (not needed for completely linear models)
254!
255 j = 1 ! counter for current nonzero
256 do i = 1, nn
257 colsta(i) = j
258 rowno(j) = 1
259 nlflag(j) = 1
260 j = j + 1
261 rowno(j) = 2
262 value(j) = 1.d0
263 nlflag(j) = 0
264 j = j + 1
265 enddo
266 colsta(nn+1) = j
267
268 qp_readmatrix = 0 ! Return value means OK
269
270end Function qp_readmatrix
271!
272!==========================================================================
273! Compute nonlinear terms and non-constant Jacobian elements
274!
275
276!> Compute nonlinear terms and non-constant Jacobian elements
277!!
278!! @include{doc} fdeval_params.dox
279Integer Function qp_fdeval( x, g, jac, rowno, jcnm, mode, ignerr, errcnt, &
280 n, nz, thread, usrmem )
281#if defined(itl)
282!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: QP_FDEval
283#endif
284 Use qp_data_sp
285 implicit none
286 integer, intent (in) :: n ! number of variables
287 integer, intent (in) :: rowno ! number of the row to be evaluated
288 integer, intent (in) :: nz ! number of nonzeros in this row
289 real*8, intent (in), dimension(n) :: x ! vector of current solution values
290 real*8, intent (in out) :: g ! constraint value
291 real*8, intent (in out), dimension(n) :: jac ! vector of derivatives for current constraint
292 integer, intent (in), dimension(nz) :: jcnm ! list of variables that appear nonlinearly
293 ! in this row. Ffor information only.
294 integer, intent (in) :: mode ! evaluation mode: 1 = function value
295 ! 2 = derivatives, 3 = both
296 integer, intent (in) :: ignerr ! if 1 then errors can be ignored as long
297 ! as errcnt is incremented
298 integer, intent (in out) :: errcnt ! error counter to be incremented in case
299 ! of function evaluation errors.
300 integer, intent (in) :: thread
301 real*8 usrmem(*) ! optional user memory
302 Integer :: i,j, k
303 real*8 :: sum
304!
305! Row 1: The objective function is nonlinear
306!
307 if ( rowno .eq. 1 ) then
308!
309! Mode = 1 or 3. Function value: G = x * Q * x / 2
310!
311 if ( mode .eq. 1 .or. mode .eq. 3 ) then
312 sum = 0.d0
313 do k = 1, nq
314 i = qrow(k)
315 j = qcol(k)
316 if ( i == j ) then
317 sum = sum + (x(i)-target(i))*q(k)*(x(j)-target(j)) / 2.d0
318 else
319 sum = sum + (x(i)-target(i))*q(k)*(x(j)-target(j))
320 endif
321 enddo
322 g = sum
323 endif
324!
325! Mode = 2 or 3: Derivative values: Q*x
326!
327 if ( mode .eq. 2 .or. mode .eq. 3 ) then
328 do i = 1, nn
329 jac(i) = 0.0d0
330 enddo
331 do k = 1, nq
332 i = qrow(k)
333 j = qcol(k)
334 if ( i == j ) then
335 jac(i) = jac(i) + q(k) * (x(i)-target(i))
336 else
337 jac(i) = jac(i) + q(k) * (x(j)-target(j))
338 jac(j) = jac(j) + q(k) * (x(i)-target(i))
339 endif
340 enddo
341 endif
342!
343! Row = 2: The row is linear and will not be called.
344!
345 endif
346 qp_fdeval = 0
347
348end Function qp_fdeval
349
350
351!> Specify the structure of the Lagrangian of the Hessian
352!!
353!! @include{doc} 2DLagrStr_params.dox
354Integer Function qp_2dlagrstr( ROWNO, COLNO, NODRV, N, M, NHESS, USRMEM )
355#if defined(itl)
356!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: QP_2DLagrStr
357#endif
358 Use qp_data_sp
359 implicit none
360 INTEGER n, m, nhess, nodrv
361 real*8 usrmem(*)
362 INTEGER rowno(nhess), colno(nhess)
363
364 Integer :: j
365
366 Do j = 1, nq
367 rowno(j) = qrow(j)
368 colno(j) = qcol(j)
369 Enddo
370 qp_2dlagrstr = 0
371
372End function qp_2dlagrstr
373
374
375!> Compute the Lagrangian of the Hessian
376!!
377!! @include{doc} 2DLagrVal_params.dox
378Integer Function qp_2dlagrval( X, U, ROWNO, COLNO, VALUE, NODRV, N, M, NHESS, USRMEM )
379#if defined(itl)
380!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: QP_2DLagrVal
381#endif
382 Use qp_data_sp
383 implicit none
384 INTEGER n, m, nhess, nodrv
385 real*8 x(n), u(m), value(nhess), usrmem(*)
386 INTEGER rowno(nhess), colno(nhess)
387
388 Integer :: j
389
390 Do j = 1, nq
391 value(j) = q(j)*u(1)
392 Enddo
393 qp_2dlagrval = 0
394
395End 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_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