CONOPT
Loading...
Searching...
No Matches
qp2.f90
Go to the documentation of this file.
1!> @file qp2.f90
2!! @ingroup FORT1THREAD_EXAMPLES
3!! @brief Similar to qp1 but uses directional 2nd derivatives
4!!
5!! For more information about the individual callbacks, please have a look at the source code.
6
7#if defined(_WIN32) && !defined(_WIN64)
8#define dec_directives_win32
9#endif
10
11!> Main program. A simple setup and call of CONOPT for a QP model
12!!
13Program qp
14
16 Use conopt
17 Use qp_data
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_2ddir ! Directional 2nd derivative routine
26 Integer, External :: std_status ! Standard callback for displaying solution status
27 Integer, External :: std_solution ! Standard callback for displaying solution values
28 Integer, External :: std_message ! Standard callback for managing messages
29 Integer, External :: std_errmsg ! Standard callback for managing error messages
30#ifdef dec_directives_win32
31!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: QP_ReadMatrix
32!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: QP_FDEval
33!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: QP_2DDir
34!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Status
35!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Solution
36!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Message
37!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_ErrMsg
38#endif
39!
40! Control vector
41!
42 INTEGER, Dimension(:), Pointer :: cntvect
43 INTEGER :: coi_error
44!
45! Local variables used to define Q
46!
47 Integer :: i,j
48
49 call startup
50!
51! Create and initialize a Control Vector
52!
53 coi_error = coi_create( cntvect )
54!
55! Tell CONOPT about the size of the model by populating the Control Vector:
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.
58!
59 j = 1
60 do i = 1, nn
61 q(j) = 1.d0
62 qrow(j) = i
63 qcol(j) = i
64 j = j + 1
65 if ( i < nn ) then
66 q(j) = 0.1d0
67 qrow(j) = i+1
68 qcol(j) = i
69 j = j + 1
70 endif
71 enddo
72 do i = 1, nn
73 target(i) = 10.d0
74 enddo
75
76 coi_error = max( coi_error, coidef_numvar( cntvect, nn ) ) ! NN variables
77 coi_error = max( coi_error, coidef_numcon( cntvect, 2 ) ) ! 1 constraint + 1 objective
78 coi_error = max( coi_error, coidef_numnz( cntvect, 2*nn ) ) ! 2*NN nonzeros in the Jacobian
79 coi_error = max( coi_error, coidef_numnlnz( cntvect, nn ) ) ! NN of which are nonlinear
80 coi_error = max( coi_error, coidef_optdir( cntvect, -1 ) ) ! Minimize
81 coi_error = max( coi_error, coidef_objcon( cntvect, 1 ) ) ! Objective is constraint 1
82!
83! Tell CONOPT about the callback routines:
84!
85 coi_error = max( coi_error, coidef_readmatrix( cntvect, qp_readmatrix ) )
86 coi_error = max( coi_error, coidef_fdeval( cntvect, qp_fdeval ) )
87 coi_error = max( coi_error, coidef_2ddir( cntvect, qp_2ddir ) )
88 coi_error = max( coi_error, coidef_status( cntvect, std_status ) )
89 coi_error = max( coi_error, coidef_solution( cntvect, std_solution ) )
90 coi_error = max( coi_error, coidef_message( cntvect, std_message ) )
91 coi_error = max( coi_error, coidef_errmsg( cntvect, std_errmsg ) )
92!
93! Default options file for debugging.
94!
95 coi_error = max( coi_error, coidef_optfile( cntvect, 'qp2.opt' ) )
96
97#if defined(CONOPT_LICENSE_INT_1) && defined(CONOPT_LICENSE_INT_2) && defined(CONOPT_LICENSE_INT_3) && defined(CONOPT_LICENSE_TEXT)
98 coi_error = max( coi_error, coidef_license( cntvect, conopt_license_int_1, conopt_license_int_2, conopt_license_int_3, conopt_license_text) )
99#endif
100
101 If ( coi_error .ne. 0 ) THEN
102 write(*,*)
103 write(*,*) '**** Fatal Error while loading CONOPT Callback routines.'
104 write(*,*)
105 call flog( "Skipping Solve due to setup errors", 1 )
106 ENDIF
107!
108! Save the solution so we can check the duals:
109!
110 do_allocate = .true.
111!
112! Start CONOPT:
113!
114 coi_error = coi_solve( cntvect )
115
116 write(*,*)
117 write(*,*) 'End of QP example 2. Return code=',coi_error
118
119 If ( coi_error /= 0 ) then
120 call flog( "Errors encountered during solution", 1 )
121 elseif ( stacalls == 0 .or. solcalls == 0 ) then
122 call flog( "Status or Solution routine was not called", 1 )
123 elseif ( sstat /= 1 .or. mstat /= 2 ) then
124 call flog( "Solver and Model Status was not as expected (1,2)", 1 )
125 Else
126 Call checkdual( 'QP 2', minimize )
127 endif
128
129 if ( coi_free(cntvect) /= 0 ) call flog( "Error while freeing control vector",1)
130
131 call flog( "Successful Solve", 0 )
133! Free solution memory
134!
135 call finalize
136
137End Program qp
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
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!> Compute nonlinear terms and non-constant Jacobian elements
228!!
229!! @include{doc} fdeval_params.dox
230Integer Function qp_fdeval( x, g, jac, rowno, jcnm, mode, ignerr, errcnt, &
231 n, nz, thread, usrmem )
232#ifdef dec_directives_win32
233!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: QP_FDEval
234#endif
235 Use qp_data
236 implicit none
237 integer, intent (in) :: n ! number of variables
238 integer, intent (in) :: rowno ! number of the row to be evaluated
239 integer, intent (in) :: nz ! number of nonzeros in this row
240 real*8, intent (in), dimension(n) :: x ! vector of current solution values
241 real*8, intent (in out) :: g ! constraint value
242 real*8, intent (in out), dimension(n) :: jac ! vector of derivatives for current constraint
243 integer, intent (in), dimension(nz) :: jcnm ! list of variables that appear nonlinearly
244 ! in this row. Ffor information only.
245 integer, intent (in) :: mode ! evaluation mode: 1 = function value
246 ! 2 = derivatives, 3 = both
247 integer, intent (in) :: ignerr ! if 1 then errors can be ignored as long
248 ! as errcnt is incremented
249 integer, intent (in out) :: errcnt ! error counter to be incremented in case
250 ! of function evaluation errors.
251 integer, intent (in) :: thread
252 real*8 usrmem(*) ! optional user memory
253 Integer :: i, j, k
254 real*8 :: sum
255!
256! Row 1: The objective function is nonlinear
257!
258 if ( rowno .eq. 1 ) then
259!
260! Mode = 1 or 3. Function value: G = x * Q * x / 2
261!
262 if ( mode .eq. 1 .or. mode .eq. 3 ) then
263 sum = 0.d0
264 do k = 1, nq
265 i = qrow(k)
266 j = qcol(k)
267 if ( i == j ) then
268 sum = sum + (x(i)-target(i))*q(k)*(x(j)-target(j))
269 else
270 sum = sum + 2*(x(i)-target(i))*q(k)*(x(j)-target(j))
271 endif
272 enddo
273 g = sum / 2.d0
274 endif
275!
276! Mode = 2 or 3: Derivative values: Q*x
277!
278 if ( mode .eq. 2 .or. mode .eq. 3 ) then
279 do i = 1, nn
280 jac(i) = 0
281 enddo
282 do k = 1, nq
283 i = qrow(k)
284 j = qcol(k)
285 if ( i == j ) then
286 jac(i) = jac(i) + q(k) * (x(i)-target(i))
287 else
288 jac(i) = jac(i) + q(k) * (x(j)-target(j))
289 jac(j) = jac(j) + q(k) * (x(i)-target(i))
290 endif
291 enddo
292 endif
293!
294! Row = 2: The row is linear and will not be called.
295!
296 endif
297 qp_fdeval = 0
298
299end Function qp_fdeval
300
301!> Computes the second derivative of a constraint in a direction
302!!
303!! @include{doc} 2DDir_params.dox
304Integer Function qp_2ddir( X, DX, D2G, ROWNO, JCNM, NODRV, N, NJ, THREAD, USRMEM )
305#ifdef dec_directives_win32
306!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: QP_2DDir
307#endif
308 Use qp_data
309 implicit none
310 INTEGER, Intent(IN) :: rowno, n, nj, thread
311 INTEGER, Intent(IN OUT) :: nodrv
312 INTEGER, Dimension(NJ), Intent(IN) :: jcnm
313 real*8, Dimension(N), Intent(IN) :: x
314 real*8, Dimension(N), Intent(IN) :: dx
315 real*8, Dimension(N), Intent(OUT) :: d2g
316 real*8, Intent(IN OUT) :: usrmem(*)
317 integer :: i,j, k
318!
319! Is only called for Irow = 1, the nonlinear objective
320! Return H*Dx = Q*Dx in D2G
321!
322 if ( rowno .eq. 1 ) then
323 do i = 1, nn
324 d2g(i) = 0.d0
325 enddo
326 do k = 1, nq
327 i = qrow(k)
328 j = qcol(k)
329 if ( i == j ) then
330 d2g(i) = d2g(i) + q(k) * dx(i)
331 else
332 d2g(i) = d2g(i) + q(k) * dx(j)
333 d2g(j) = d2g(j) + q(k) * dx(i)
334 endif
335 enddo
336 qp_2ddir = 0
337 else
338 qp_2ddir = 1
339 endif
340
341End Function qp_2ddir
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
subroutine checkdual(case, minmax)
Definition comdecl.f90:432
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_optfile(cntvect, optfile)
define callback routine for defining an options file.
Definition conopt.f90:928
integer(c_int) function coidef_2ddir(cntvect, coi_2ddir)
define callback routine for computing the second derivative for a constraint in a direction.
Definition conopt.f90:1421
integer(c_int) function coidef_license(cntvect, licint1, licint2, licint3, licstring)
define the License Information.
Definition conopt.f90:293
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_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
#define nj
Definition mp_trans.c:46
integer solcalls
Definition comdecl.f90:15
integer sstat
Definition comdecl.f90:18
subroutine finalize
Definition comdecl.f90:79
integer, parameter minimize
Definition comdecl.f90:31
integer stacalls
Definition comdecl.f90:14
subroutine flog(msg, code)
Definition comdecl.f90:62
logical do_allocate
Definition comdecl.f90:27
integer mstat
Definition comdecl.f90:17
subroutine startup
Definition comdecl.f90:41
real(8) target
integer, parameter nn
real(8), dimension(nn, nn) q
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_2ddir(x, dx, d2g, rowno, jcnm, nodrv, n, nj, thread, usrmem)
Computes the second derivative of a constraint in a direction.
Definition qp2.f90:289