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#if defined(_WIN32) && !defined(_WIN64)
11#define dec_directives_win32
12#endif
13
14!> Main program. A simple setup and call of CONOPT for a QP model
15!!
16Program qp
17
19 Use conopt
20 Use qp_data_sp
21 implicit None
22!
23! Declare the user callback routines as Integer, External:
24!
25 Integer, External :: qp_readmatrix ! Mandatory Matrix definition routine defined below
26 Integer, External :: qp_fdeval ! Function and Derivative evaluation routine
27 ! needed a nonlinear model.
28 Integer, External :: qp_2dlagrstr ! 2nd derivative routine for structure
29 Integer, External :: qp_2dlagrval ! 2nd derivative routine for values
30 Integer, External :: std_status ! Standard callback for displaying solution status
31 Integer, External :: std_solution ! Standard callback for displaying solution values
32 Integer, External :: std_message ! Standard callback for managing messages
33 Integer, External :: std_errmsg ! Standard callback for managing error messages
34#ifdef dec_directives_win32
35!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: QP_ReadMatrix
36!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: QP_FDEval
37!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: QP_2DLagrStr
38!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: QP_2DLagrVal
39!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Status
40!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Solution
41!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Message
42!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_ErrMsg
43#endif
44!
45! Control vector
46!
47 INTEGER, Dimension(:), Pointer :: cntvect
48 INTEGER :: coi_error
49!
50! Local variables used to define Q
51!
52 Integer :: i,j,k
53
54 call startup
55!
56! We will create a QP model with one constraint, sum(i, x(i) ) = 1
57! and NN variables. NN and the Q matrix are defined in Module QP_Data_Sp.
58!
59 j = 0
60 do i = 1, nn
61 j = j + 1
62 q(j) = 1.d0
63 qrow(j) = i
64 qcol(j) = i
65 do k = i+1, nn, 2
66 j = j + 1
67 q(j) = 0.5d0/nn ! Make sure the matrix is diagonally dominant.
68 qrow(j) = k
69 qcol(j) = i
70 enddo
71 enddo
72 nq = j
73 do i = 1, nn
74 target(i) = 0.1d0
75 enddo
76!
77! Create and initialize a Control Vector
78!
79 coi_error = coi_create( 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(CONOPT_LICENSE_INT_1) && defined(CONOPT_LICENSE_INT_2) && defined(CONOPT_LICENSE_INT_3) && defined(CONOPT_LICENSE_TEXT)
105 coi_error = max( coi_error, coidef_license( cntvect, conopt_license_int_1, conopt_license_int_2, conopt_license_int_3, conopt_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!
178! Free solution memory
179!
181
182End Program qp
183!
184! ============================================================================
185! Define information about the model:
186!
187
188!> Define information about the model
189!!
190!! @include{doc} readMatrix_params.dox
191Integer Function qp_readmatrix( lower, curr, upper, vsta, type, rhs, esta, &
192 colsta, rowno, value, nlflag, n, m, nz, &
193 usrmem )
194#ifdef dec_directives_win32
195!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: QP_ReadMatrix
196#endif
197 Use qp_data_sp
198 implicit none
199 integer, intent (in) :: n ! number of variables
200 integer, intent (in) :: m ! number of constraints
201 integer, intent (in) :: nz ! number of nonzeros
202 real*8, intent (in out), dimension(n) :: lower ! vector of lower bounds
203 real*8, intent (in out), dimension(n) :: curr ! vector of initial values
204 real*8, intent (in out), dimension(n) :: upper ! vector of upper bounds
205 integer, intent (in out), dimension(n) :: vsta ! vector of initial variable status
206 ! (not defined here)
207 integer, intent (out), dimension(m) :: type ! vector of equation types
208 integer, intent (in out), dimension(m) :: esta ! vector of initial equation status
209 ! (not defined here)
210 real*8, intent (in out), dimension(m) :: rhs ! vector of right hand sides
211 integer, intent (in out), dimension(n+1) :: colsta ! vector with start of column indices
212 integer, intent (out), dimension(nz) :: rowno ! vector of row numbers
213 integer, intent (in out), dimension(nz) :: nlflag ! vector of nonlinearity flags
214 real*8, intent (in out), dimension(nz) :: value ! vector of matrix values
215 real*8 usrmem(*) ! optional user memory
216 integer :: i, j
217!
218! Information about Variables:
219! Default: Lower = -Inf, Curr = 0, and Upper = +inf.
220! Default: the status information in Vsta is not used.
221!
222! Nondefault: Lower bound = 0 for all variables
223!
224 do i = 1, nn
225 lower(i) = 0.0d0
226 enddo
227!
228! Information about Constraints:
229! Default: Rhs = 0
230! Default: the status information in Esta and the function
231! value in FV are not used.
232! Default: Type: There is no default.
233! 0 = Equality,
234! 1 = Greater than or equal,
235! 2 = Less than or equal,
236! 3 = Non binding.
237!
238! Constraint 1 (Objective)
239! Rhs = 0.0 and type Non binding
240!
241 type(1) = 3
242!
243! Constraint 2 (Sum of all vars = 1)
244! Rhs = 1 and type Equality
245!
246 rhs(2) = 1.d0
247 type(2) = 0
248!
249! Information about the Jacobian. CONOPT expects a columnwise
250! representation in Rowno, Value, Nlflag and Colsta.
251!
252! Colsta = Start of column indices (No Defaults):
253! Rowno = Row indices
254! Value = Value of derivative (by default only linear
255! derivatives are used)
256! Nlflag = 0 for linear and 1 for nonlinear derivative
257! (not needed for completely linear models)
258!
259 j = 1 ! counter for current nonzero
260 do i = 1, nn
261 colsta(i) = j
262 rowno(j) = 1
263 nlflag(j) = 1
264 j = j + 1
265 rowno(j) = 2
266 value(j) = 1.d0
267 nlflag(j) = 0
268 j = j + 1
269 enddo
270 colsta(nn+1) = j
271
272 qp_readmatrix = 0 ! Return value means OK
273
274end Function qp_readmatrix
275!
276!==========================================================================
277! Compute nonlinear terms and non-constant Jacobian elements
278!
279
280!> Compute nonlinear terms and non-constant Jacobian elements
281!!
282!! @include{doc} fdeval_params.dox
283Integer Function qp_fdeval( x, g, jac, rowno, jcnm, mode, ignerr, errcnt, &
284 n, nz, thread, usrmem )
285#ifdef dec_directives_win32
286!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: QP_FDEval
287#endif
288 Use qp_data_sp
289 implicit none
290 integer, intent (in) :: n ! number of variables
291 integer, intent (in) :: rowno ! number of the row to be evaluated
292 integer, intent (in) :: nz ! number of nonzeros in this row
293 real*8, intent (in), dimension(n) :: x ! vector of current solution values
294 real*8, intent (in out) :: g ! constraint value
295 real*8, intent (in out), dimension(n) :: jac ! vector of derivatives for current constraint
296 integer, intent (in), dimension(nz) :: jcnm ! list of variables that appear nonlinearly
297 ! in this row. Ffor information only.
298 integer, intent (in) :: mode ! evaluation mode: 1 = function value
299 ! 2 = derivatives, 3 = both
300 integer, intent (in) :: ignerr ! if 1 then errors can be ignored as long
301 ! as errcnt is incremented
302 integer, intent (in out) :: errcnt ! error counter to be incremented in case
303 ! of function evaluation errors.
304 integer, intent (in) :: thread
305 real*8 usrmem(*) ! optional user memory
306 Integer :: i,j, k
307 real*8 :: sum
308!
309! Row 1: The objective function is nonlinear
310!
311 if ( rowno .eq. 1 ) then
312!
313! Mode = 1 or 3. Function value: G = x * Q * x / 2
314!
315 if ( mode .eq. 1 .or. mode .eq. 3 ) then
316 sum = 0.d0
317 do k = 1, nq
318 i = qrow(k)
319 j = qcol(k)
320 if ( i == j ) then
321 sum = sum + (x(i)-target(i))*q(k)*(x(j)-target(j)) / 2.d0
322 else
323 sum = sum + (x(i)-target(i))*q(k)*(x(j)-target(j))
324 endif
325 enddo
326 g = sum
327 endif
328!
329! Mode = 2 or 3: Derivative values: Q*x
330!
331 if ( mode .eq. 2 .or. mode .eq. 3 ) then
332 do i = 1, nn
333 jac(i) = 0.0d0
334 enddo
335 do k = 1, nq
336 i = qrow(k)
337 j = qcol(k)
338 if ( i == j ) then
339 jac(i) = jac(i) + q(k) * (x(i)-target(i))
340 else
341 jac(i) = jac(i) + q(k) * (x(j)-target(j))
342 jac(j) = jac(j) + q(k) * (x(i)-target(i))
343 endif
344 enddo
345 endif
346!
347! Row = 2: The row is linear and will not be called.
348!
349 endif
350 qp_fdeval = 0
351
352end Function qp_fdeval
353
354
355!> Specify the structure of the Lagrangian of the Hessian
356!!
357!! @include{doc} 2DLagrStr_params.dox
358Integer Function qp_2dlagrstr( ROWNO, COLNO, NODRV, N, M, NHESS, USRMEM )
359#ifdef dec_directives_win32
360!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: QP_2DLagrStr
361#endif
362 Use qp_data_sp
363 implicit none
364 INTEGER n, m, nhess, nodrv
365 real*8 usrmem(*)
366 INTEGER rowno(nhess), colno(nhess)
367
368 Integer :: j
369
370 Do j = 1, nq
371 rowno(j) = qrow(j)
372 colno(j) = qcol(j)
373 Enddo
374 qp_2dlagrstr = 0
375
376End function qp_2dlagrstr
377
378
379!> Compute the Lagrangian of the Hessian
380!!
381!! @include{doc} 2DLagrVal_params.dox
382Integer Function qp_2dlagrval( X, U, ROWNO, COLNO, VALUE, NODRV, N, M, NHESS, USRMEM )
383#ifdef dec_directives_win32
384!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: QP_2DLagrVal
385#endif
386 Use qp_data_sp
387 implicit none
388 INTEGER n, m, nhess, nodrv
389 real*8 x(n), u(m), value(nhess), usrmem(*)
390 INTEGER rowno(nhess), colno(nhess)
391
392 Integer :: j
393
394 Do j = 1, nq
395 value(j) = q(j)*u(1)
396 Enddo
397 qp_2dlagrval = 0
398
399End function qp_2dlagrval
integer function std_solution(xval, xmar, xbas, xsta, yval, ymar, ybas, ysta, n, m, usrmem)
Definition comdecl.f90:170
integer function std_status(modsta, solsta, iter, objval, usrmem)
Definition comdecl.f90:126
integer function std_message(smsg, dmsg, nmsg, llen, usrmem, msgv)
Definition comdecl.f90:243
integer function std_errmsg(rowno, colno, posno, msglen, usrmem, msg)
Definition comdecl.f90:286
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_optfile(cntvect, optfile)
define callback routine for defining an options file.
Definition conopt.f90:928
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_hessfac(cntvect, hessfac)
factor for Hessian density relative to Jacobian density HessFac.
Definition conopt.f90:858
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
integer solcalls
Definition comdecl.f90:15
integer sstat
Definition comdecl.f90:18
subroutine finalize
Definition comdecl.f90:79
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
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:204
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:288
integer function qp_2dlagrval(x, u, rowno, colno, value, nodrv, n, m, nhess, usrmem)
Compute the Lagrangian of the Hessian.
Definition qp3.f90:302
integer function qp_2dlagrstr(rowno, colno, nodrv, n, m, nhess, usrmem)
Specify the structure of the Lagrangian of the Hessian.
Definition qp3.f90:282