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