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