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!> Main program. A simple setup and call of CONOPT for a QP model
8!!
9Program qp
10
11 Use proginfo
12 Use coidef
13 Use qp_data
14 implicit None
15!
16! Declare the user callback routines as Integer, External:
17!
18 Integer, External :: qp_readmatrix ! Mandatory Matrix definition routine defined below
19 Integer, External :: qp_fdeval ! Function and Derivative evaluation routine
20 ! needed a nonlinear model.
21 Integer, External :: qp_2ddir ! Directional 2nd derivative routine
22 Integer, External :: std_status ! Standard callback for displaying solution status
23 Integer, External :: std_solution ! Standard callback for displaying solution values
24 Integer, External :: std_message ! Standard callback for managing messages
25 Integer, External :: std_errmsg ! Standard callback for managing error messages
26#if defined(itl)
27!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: QP_ReadMatrix
28!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: QP_FDEval
29!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: QP_2DDir
30!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Status
31!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Solution
32!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Message
33!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_ErrMsg
34#endif
35!
36! Control vector
37!
38 INTEGER :: numcallback
39 INTEGER, Dimension(:), Pointer :: cntvect
40 INTEGER :: coi_error
41!
42! Local variables used to define Q
43!
44 Integer :: i,j
45
46 call startup
47!
48! Create and initialize a Control Vector
49!
50 numcallback = coidef_size()
51 Allocate( cntvect(numcallback) )
52 coi_error = coidef_inifort( cntvect )
53!
54! Tell CONOPT about the size of the model by populating the Control Vector:
55! We will create a QP model with one constraint, sum(i, x(i) ) = 1
56! and NN variables. NN and the Q matrix are defined in Module QP_Data.
57!
58 j = 1
59 do i = 1, nn
60 q(j) = 1.d0
61 qrow(j) = i
62 qcol(j) = i
63 j = j + 1
64 if ( i < nn ) then
65 q(j) = 0.1d0
66 qrow(j) = i+1
67 qcol(j) = i
68 j = j + 1
69 endif
70 enddo
71 do i = 1, nn
72 target(i) = 10.d0
73 enddo
74
75 coi_error = max( coi_error, coidef_numvar( cntvect, nn ) ) ! NN variables
76 coi_error = max( coi_error, coidef_numcon( cntvect, 2 ) ) ! 1 constraint + 1 objective
77 coi_error = max( coi_error, coidef_numnz( cntvect, 2*nn ) ) ! 2*NN nonzeros in the Jacobian
78 coi_error = max( coi_error, coidef_numnlnz( cntvect, nn ) ) ! NN of which are nonlinear
79 coi_error = max( coi_error, coidef_optdir( cntvect, -1 ) ) ! Minimize
80 coi_error = max( coi_error, coidef_objcon( cntvect, 1 ) ) ! Objective is constraint 1
81!
82! Tell CONOPT about the callback routines:
83!
84 coi_error = max( coi_error, coidef_readmatrix( cntvect, qp_readmatrix ) )
85 coi_error = max( coi_error, coidef_fdeval( cntvect, qp_fdeval ) )
86 coi_error = max( coi_error, coidef_2ddir( cntvect, qp_2ddir ) )
87 coi_error = max( coi_error, coidef_status( cntvect, std_status ) )
88 coi_error = max( coi_error, coidef_solution( cntvect, std_solution ) )
89 coi_error = max( coi_error, coidef_message( cntvect, std_message ) )
90 coi_error = max( coi_error, coidef_errmsg( cntvect, std_errmsg ) )
91!
92! Default options file for debugging.
93!
94 coi_error = max( coi_error, coidef_optfile( cntvect, 'qp2.opt' ) )
95
96#if defined(LICENSE_INT_1) && defined(LICENSE_INT_2) && defined(LICENSE_INT_3) && defined(LICENSE_TEXT)
97 coi_error = max( coi_error, coidef_license( cntvect, license_int_1, license_int_2, license_int_3, license_text) )
98#endif
99
100 If ( coi_error .ne. 0 ) THEN
101 write(*,*)
102 write(*,*) '**** Fatal Error while loading CONOPT Callback routines.'
103 write(*,*)
104 call flog( "Skipping Solve due to setup errors", 1 )
105 ENDIF
106!
107! Save the solution so we can check the duals:
108!
109 do_allocate = .true.
110!
111! Start CONOPT:
112!
113 coi_error = coi_solve( cntvect )
114
115 write(*,*)
116 write(*,*) 'End of QP example 2. Return code=',coi_error
117
118 If ( coi_error /= 0 ) then
119 call flog( "Errors encountered during solution", 1 )
120 elseif ( stacalls == 0 .or. solcalls == 0 ) then
121 call flog( "Status or Solution routine was not called", 1 )
122 elseif ( sstat /= 1 .or. mstat /= 2 ) then
123 call flog( "Solver and Model Status was not as expected (1,2)", 1 )
124 Else
125 Call checkdual( 'QP 2', minimize )
126 endif
127
128 if ( coi_free(cntvect) /= 0 ) call flog( "Error while freeing control vector",1)
129
130 call flog( "Successful Solve", 0 )
131
132End Program qp
133
134!> Define information about the model
135!!
136!! @include{doc} readMatrix_params.dox
137Integer Function qp_readmatrix( lower, curr, upper, vsta, type, rhs, esta, &
138 colsta, rowno, value, nlflag, n, m, nz, &
139 usrmem )
140#if defined(itl)
141!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: QP_ReadMatrix
142#endif
143 Use qp_data
144 implicit none
145 integer, intent (in) :: n ! number of variables
146 integer, intent (in) :: m ! number of constraints
147 integer, intent (in) :: nz ! number of nonzeros
148 real*8, intent (in out), dimension(n) :: lower ! vector of lower bounds
149 real*8, intent (in out), dimension(n) :: curr ! vector of initial values
150 real*8, intent (in out), dimension(n) :: upper ! vector of upper bounds
151 integer, intent (in out), dimension(n) :: vsta ! vector of initial variable status
152 ! (not defined here)
153 integer, intent (out), dimension(m) :: type ! vector of equation types
154 integer, intent (in out), dimension(m) :: esta ! vector of initial equation status
155 ! (not defined here)
156 real*8, intent (in out), dimension(m) :: rhs ! vector of right hand sides
157 integer, intent (in out), dimension(n+1) :: colsta ! vector with start of column indices
158 integer, intent (out), dimension(nz) :: rowno ! vector of row numbers
159 integer, intent (in out), dimension(nz) :: nlflag ! vector of nonlinearity flags
160 real*8, intent (in out), dimension(nz) :: value ! vector of matrix values
161 real*8 usrmem(*) ! optional user memory
162 integer :: i, j
163!
164! Information about Variables:
165! Default: Lower = -Inf, Curr = 0, and Upper = +inf.
166! Default: the status information in Vsta is not used.
167!
168! Nondefault: Lower bound = 0 for all variables
169!
170 do i = 1, nn
171 lower(i) = 0.0d0
172 enddo
173!
174! Information about Constraints:
175! Default: Rhs = 0
176! Default: the status information in Esta and the function
177! value in FV are not used.
178! Default: Type: There is no default.
179! 0 = Equality,
180! 1 = Greater than or equal,
181! 2 = Less than or equal,
182! 3 = Non binding.
183!
184! Constraint 1 (Objective)
185! Rhs = 0.0 and type Non binding
186!
187 type(1) = 3
188!
189! Constraint 2 (Sum of all vars = 1)
190! Rhs = 1 and type Equality
191!
192 rhs(2) = 1.d0
193 type(2) = 0
194!
195! Information about the Jacobian. We use the standard method with
196! Rowno, Value, Nlflag and Colsta and we do not use Colno.
197!
198! Colsta = Start of column indices (No Defaults):
199! Rowno = Row indices
200! Value = Value of derivative (by default only linear
201! derivatives are used)
202! Nlflag = 0 for linear and 1 for nonlinear derivative
203! (not needed for completely linear models)
204!
205 j = 1 ! counter for current nonzero
206 do i = 1, nn
207 colsta(i) = j
208 rowno(j) = 1
209 nlflag(j) = 1
210 j = j + 1
211 rowno(j) = 2
212 value(j) = 1.d0
213 nlflag(j) = 0
214 j = j + 1
215 enddo
216 colsta(nn+1) = j
217
218 qp_readmatrix = 0 ! Return value means OK
219
220end Function qp_readmatrix
221
222!> Compute nonlinear terms and non-constant Jacobian elements
223!!
224!! @include{doc} fdeval_params.dox
225Integer Function qp_fdeval( x, g, jac, rowno, jcnm, mode, ignerr, errcnt, &
226 n, nz, thread, usrmem )
227#if defined(itl)
228!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: QP_FDEval
229#endif
230 Use qp_data
231 implicit none
232 integer, intent (in) :: n ! number of variables
233 integer, intent (in) :: rowno ! number of the row to be evaluated
234 integer, intent (in) :: nz ! number of nonzeros in this row
235 real*8, intent (in), dimension(n) :: x ! vector of current solution values
236 real*8, intent (in out) :: g ! constraint value
237 real*8, intent (in out), dimension(n) :: jac ! vector of derivatives for current constraint
238 integer, intent (in), dimension(nz) :: jcnm ! list of variables that appear nonlinearly
239 ! in this row. Ffor information only.
240 integer, intent (in) :: mode ! evaluation mode: 1 = function value
241 ! 2 = derivatives, 3 = both
242 integer, intent (in) :: ignerr ! if 1 then errors can be ignored as long
243 ! as errcnt is incremented
244 integer, intent (in out) :: errcnt ! error counter to be incremented in case
245 ! of function evaluation errors.
246 integer, intent (in) :: thread
247 real*8 usrmem(*) ! optional user memory
248 Integer :: i, j, k
249 real*8 :: sum
250!
251! Row 1: The objective function is nonlinear
252!
253 if ( rowno .eq. 1 ) then
254!
255! Mode = 1 or 3. Function value: G = x * Q * x / 2
256!
257 if ( mode .eq. 1 .or. mode .eq. 3 ) then
258 sum = 0.d0
259 do k = 1, nq
260 i = qrow(k)
261 j = qcol(k)
262 if ( i == j ) then
263 sum = sum + (x(i)-target(i))*q(k)*(x(j)-target(j))
264 else
265 sum = sum + 2*(x(i)-target(i))*q(k)*(x(j)-target(j))
266 endif
267 enddo
268 g = sum / 2.d0
269 endif
270!
271! Mode = 2 or 3: Derivative values: Q*x
272!
273 if ( mode .eq. 2 .or. mode .eq. 3 ) then
274 do i = 1, nn
275 jac(i) = 0
276 enddo
277 do k = 1, nq
278 i = qrow(k)
279 j = qcol(k)
280 if ( i == j ) then
281 jac(i) = jac(i) + q(k) * (x(i)-target(i))
282 else
283 jac(i) = jac(i) + q(k) * (x(j)-target(j))
284 jac(j) = jac(j) + q(k) * (x(i)-target(i))
285 endif
286 enddo
287 endif
288!
289! Row = 2: The row is linear and will not be called.
290!
291 endif
292 qp_fdeval = 0
293
294end Function qp_fdeval
295
296!> Computes the second derivative of a constraint in a direction
297!!
298!! @include{doc} 2DDir_params.dox
299Integer Function qp_2ddir( X, DX, D2G, ROWNO, JCNM, NODRV, N, NJ, THREAD, USRMEM )
300#if defined(itl)
301!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: QP_2DDir
302#endif
303 Use qp_data
304 implicit none
305 INTEGER, Intent(IN) :: rowno, n, nj, thread
306 INTEGER, Intent(IN OUT) :: nodrv
307 INTEGER, Dimension(NJ), Intent(IN) :: jcnm
308 real*8, Dimension(N), Intent(IN) :: x
309 real*8, Dimension(N), Intent(IN) :: dx
310 real*8, Dimension(N), Intent(OUT) :: d2g
311 real*8, Intent(IN OUT) :: usrmem(*)
312 integer :: i,j, k
313!
314! Is only called for Irow = 1, the nonlinear objective
315! Return H*Dx = Q*Dx in D2G
316!
317 if ( rowno .eq. 1 ) then
318 do i = 1, nn
319 d2g(i) = 0.d0
320 enddo
321 do k = 1, nq
322 i = qrow(k)
323 j = qcol(k)
324 if ( i == j ) then
325 d2g(i) = d2g(i) + q(k) * dx(i)
326 else
327 d2g(i) = d2g(i) + q(k) * dx(j)
328 d2g(j) = d2g(j) + q(k) * dx(i)
329 endif
330 enddo
331 qp_2ddir = 0
332 else
333 qp_2ddir = 1
334 endif
335
336End Function qp_2ddir
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
subroutine checkdual(case, minmax)
Definition comdecl.f90:365
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_2ddir(cntvect, coi_2ddir)
define callback routine for computing the second derivative for a constraint in a direction.
integer function coidef_optfile(cntvect, optfile)
define callback routine for defining an options file.
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_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
#define nj
Definition mp_trans.c:46
integer solcalls
Definition comdecl.f90:9
integer sstat
Definition comdecl.f90:12
integer, parameter minimize
Definition comdecl.f90:25
integer stacalls
Definition comdecl.f90:8
subroutine flog(msg, code)
Definition comdecl.f90:56
logical do_allocate
Definition comdecl.f90:21
integer mstat
Definition comdecl.f90:11
subroutine startup
Definition comdecl.f90:35
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: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_2ddir(x, dx, d2g, rowno, jcnm, nodrv, n, nj, thread, usrmem)
Computes the second derivative of a constraint in a direction.
Definition qp2.f90:300