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