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
133End Program qp
134!
135! ============================================================================
136! Define information about the model:
137!
138
139!> Define information about the model
140!!
141!! @include{doc} readMatrix_params.dox
142Integer Function qp_readmatrix( lower, curr, upper, vsta, type, rhs, esta, &
143 colsta, rowno, value, nlflag, n, m, nz, &
144 usrmem )
145#ifdef dec_directives_win32
146!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: QP_ReadMatrix
147#endif
148 Use qp_data_sp
149 implicit none
150 integer, intent (in) :: n ! number of variables
151 integer, intent (in) :: m ! number of constraints
152 integer, intent (in) :: nz ! number of nonzeros
153 real*8, intent (in out), dimension(n) :: lower ! vector of lower bounds
154 real*8, intent (in out), dimension(n) :: curr ! vector of initial values
155 real*8, intent (in out), dimension(n) :: upper ! vector of upper bounds
156 integer, intent (in out), dimension(n) :: vsta ! vector of initial variable status
157 ! (not defined here)
158 integer, intent (out), dimension(m) :: type ! vector of equation types
159 integer, intent (in out), dimension(m) :: esta ! vector of initial equation status
160 ! (not defined here)
161 real*8, intent (in out), dimension(m) :: rhs ! vector of right hand sides
162 integer, intent (in out), dimension(n+1) :: colsta ! vector with start of column indices
163 integer, intent (out), dimension(nz) :: rowno ! vector of row numbers
164 integer, intent (in out), dimension(nz) :: nlflag ! vector of nonlinearity flags
165 real*8, intent (in out), dimension(nz) :: value ! vector of matrix values
166 real*8 usrmem(*) ! optional user memory
167 integer :: i, j
168!
169! Information about Variables:
170! Default: Lower = -Inf, Curr = 0, and Upper = +inf.
171! Default: the status information in Vsta is not used.
172!
173! Nondefault: Lower bound = 0 for all variables
174!
175 do i = 1, nn
176 lower(i) = 0.0d0
177 enddo
178!
179! Information about Constraints:
180! Default: Rhs = 0
181! Default: the status information in Esta and the function
182! value in FV are not used.
183! Default: Type: There is no default.
184! 0 = Equality,
185! 1 = Greater than or equal,
186! 2 = Less than or equal,
187! 3 = Non binding.
188!
189! Constraint 1 (Objective)
190! Rhs = 0.0 and type Non binding
191!
192 type(1) = 3
193!
194! Constraint 2 (Sum of all vars = 1)
195! Rhs = 1 and type Equality
196!
197 rhs(2) = 1.d0
198 type(2) = 0
199!
200! Information about the Jacobian. CONOPT expects a columnwise
201! representation in Rowno, Value, Nlflag and Colsta.
202!
203! Colsta = Start of column indices (No Defaults):
204! Rowno = Row indices
205! Value = Value of derivative (by default only linear
206! derivatives are used)
207! Nlflag = 0 for linear and 1 for nonlinear derivative
208! (not needed for completely linear models)
209!
210 j = 1 ! counter for current nonzero
211 do i = 1, nn
212 colsta(i) = j
213 rowno(j) = 1
214 nlflag(j) = 1
215 j = j + 1
216 rowno(j) = 2
217 value(j) = 1.d0
218 nlflag(j) = 0
219 j = j + 1
220 enddo
221 colsta(nn+1) = j
222
223 qp_readmatrix = 0 ! Return value means OK
224
225end Function qp_readmatrix
226!
227!==========================================================================
228! Compute nonlinear terms and non-constant Jacobian elements
229!
230
231!> Compute nonlinear terms and non-constant Jacobian elements
232!!
233!! @include{doc} fdeval_params.dox
234Integer Function qp_fdeval( x, g, jac, rowno, jcnm, mode, ignerr, errcnt, &
235 n, nz, thread, usrmem )
236#ifdef dec_directives_win32
237!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: QP_FDEval
238#endif
239 Use qp_data_sp
240 implicit none
241 integer, intent (in) :: n ! number of variables
242 integer, intent (in) :: rowno ! number of the row to be evaluated
243 integer, intent (in) :: nz ! number of nonzeros in this row
244 real*8, intent (in), dimension(n) :: x ! vector of current solution values
245 real*8, intent (in out) :: g ! constraint value
246 real*8, intent (in out), dimension(n) :: jac ! vector of derivatives for current constraint
247 integer, intent (in), dimension(nz) :: jcnm ! list of variables that appear nonlinearly
248 ! in this row. Ffor information only.
249 integer, intent (in) :: mode ! evaluation mode: 1 = function value
250 ! 2 = derivatives, 3 = both
251 integer, intent (in) :: ignerr ! if 1 then errors can be ignored as long
252 ! as errcnt is incremented
253 integer, intent (in out) :: errcnt ! error counter to be incremented in case
254 ! of function evaluation errors.
255 integer, intent (in) :: thread
256 real*8 usrmem(*) ! optional user memory
257 Integer :: i,j, k
258 real*8 :: sum
259!
260! Row 1: The objective function is nonlinear
261!
262 if ( rowno .eq. 1 ) then
263!
264! Mode = 1 or 3. Function value: G = x * Q * x / 2
265!
266 if ( mode .eq. 1 .or. mode .eq. 3 ) then
267 sum = 0.d0
268 do k = 1, nq
269 i = qrow(k)
270 j = qcol(k)
271 if ( i == j ) then
272 sum = sum + (x(i)-target(i))*q(k)*(x(j)-target(j)) / 2.d0
273 else
274 sum = sum + (x(i)-target(i))*q(k)*(x(j)-target(j))
275 endif
276 enddo
277 g = sum
278 endif
279!
280! Mode = 2 or 3: Derivative values: Q*x
281!
282 if ( mode .eq. 2 .or. mode .eq. 3 ) then
283 do i = 1, nn
284 jac(i) = 0.0d0
285 enddo
286 do k = 1, nq
287 i = qrow(k)
288 j = qcol(k)
289 if ( i == j ) then
290 jac(i) = jac(i) + q(k) * (x(i)-target(i))
291 else
292 jac(i) = jac(i) + q(k) * (x(j)-target(j))
293 jac(j) = jac(j) + q(k) * (x(i)-target(i))
294 endif
295 enddo
296 endif
297!
298! Row = 2: The row is linear and will not be called.
299!
300 endif
301 qp_fdeval = 0
302
303end Function qp_fdeval
304
305
306!> Specify the structure of the Lagrangian of the Hessian
307!!
308!! @include{doc} 2DLagrStr_params.dox
309Integer Function qp_2dlagrstr( ROWNO, COLNO, NODRV, N, M, NHESS, USRMEM )
310#ifdef dec_directives_win32
311!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: QP_2DLagrStr
312#endif
313 Use qp_data_sp
314 implicit none
315 INTEGER n, m, nhess, nodrv
316 real*8 usrmem(*)
317 INTEGER rowno(nhess), colno(nhess)
318
319 Integer :: j
320
321 Do j = 1, nq
322 rowno(j) = qrow(j)
323 colno(j) = qcol(j)
324 Enddo
325 qp_2dlagrstr = 0
326
327End function qp_2dlagrstr
328
329
330!> Compute the Lagrangian of the Hessian
331!!
332!! @include{doc} 2DLagrVal_params.dox
333Integer Function qp_2dlagrval( X, U, ROWNO, COLNO, VALUE, NODRV, N, M, NHESS, USRMEM )
334#ifdef dec_directives_win32
335!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: QP_2DLagrVal
336#endif
337 Use qp_data_sp
338 implicit none
339 INTEGER n, m, nhess, nodrv
340 real*8 x(n), u(m), value(nhess), usrmem(*)
341 INTEGER rowno(nhess), colno(nhess)
342
343 Integer :: j
344
345 Do j = 1, nq
346 value(j) = q(j)*u(1)
347 Enddo
348 qp_2dlagrval = 0
349
350End function qp_2dlagrval
integer function std_solution(xval, xmar, xbas, xsta, yval, ymar, ybas, ysta, n, m, usrmem)
Definition comdecl.f90:132
integer function std_status(modsta, solsta, iter, objval, usrmem)
Definition comdecl.f90:88
integer function std_message(smsg, dmsg, nmsg, llen, usrmem, msgv)
Definition comdecl.f90:205
integer function std_errmsg(rowno, colno, posno, msglen, usrmem, msg)
Definition comdecl.f90:248
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
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:200
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:284
integer function qp_2dlagrval(x, u, rowno, colno, value, nodrv, n, m, nhess, usrmem)
Compute the Lagrangian of the Hessian.
Definition qp3.f90:298
integer function qp_2dlagrstr(rowno, colno, nodrv, n, m, nhess, usrmem)
Specify the structure of the Lagrangian of the Hessian.
Definition qp3.f90:278