CONOPT
Loading...
Searching...
No Matches
qp7.f90
Go to the documentation of this file.
1!> @file qp7.f90
2!! @ingroup FORT1THREAD_EXAMPLES
3!! @brief Similar to qp5 but with the objective defined as a positive
4!! variable, i.e. it cannot be removed in the post-triangle
5!!
6!! For more information about the individual callbacks, please have a look at the source code.
7
8#if defined(_WIN32) && !defined(_WIN64)
9#define dec_directives_win32
10#endif
11
12!> Main program. A simple setup and call of CONOPT for a QP model
13!!
14Program qp
15
17 Use conopt
18 Use qp_data
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 :: 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 :: Std_Status
34!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Solution
35!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Message
36!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_ErrMsg
37#endif
38!
39! Control vector
40!
41 INTEGER, Dimension(:), Pointer :: cntvect
42 INTEGER :: coi_error
43!
44! Local variables used to define Q
45!
46 Integer :: i,j
47
48 Call startup
49!
50! Create and initialize a Control Vector
51!
52 coi_error = coi_create( 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 declared 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+1 ) ) ! NN+1 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+1 ) ) ! 2*NN+1 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_objvar( cntvect, nn+1 ) ) ! Objective is variable NN+1
81 coi_error = max( coi_error, coidef_optfile( cntvect, 'test.opt' ) )
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_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#if defined(CONOPT_LICENSE_INT_1) && defined(CONOPT_LICENSE_INT_2) && defined(CONOPT_LICENSE_INT_3) && defined(CONOPT_LICENSE_TEXT)
93 coi_error = max( coi_error, coidef_license( cntvect, conopt_license_int_1, conopt_license_int_2, conopt_license_int_3, conopt_license_text) )
94#endif
95
96 If ( coi_error .ne. 0 ) THEN
97 write(*,*)
98 write(*,*) '**** Fatal Error while loading CONOPT Callback routines.'
99 write(*,*)
100 call flog( "Skipping Solve due to setup errors", 1 )
101 ENDIF
102!
103! Start CONOPT:
104!
105 coi_error = coi_solve( cntvect )
106
107 write(*,*)
108 write(*,*) 'End of QP example 7. Return code=',coi_error
109
110 If ( coi_error /= 0 ) then
111 call flog( "Errors encountered during solution", 1 )
112 elseif ( stacalls == 0 .or. solcalls == 0 ) then
113 call flog( "Status or Solution routine was not called", 1 )
114 elseif ( sstat /= 1 .or. mstat /= 2 ) then
115 call flog( "Solver and Model Status was not as expected (1,2)", 1 )
116 endif
117
118 if ( coi_free(cntvect) /= 0 ) call flog( "Error while freeing control vector",1)
119
120 call flog( "Successful Solve", 0 )
121
122End Program qp
123
124!> Define information about the model
125!!
126!! @include{doc} readMatrix_params.dox
127Integer Function qp_readmatrix( lower, curr, upper, vsta, type, rhs, esta, &
128 colsta, rowno, value, nlflag, n, m, nz, &
129 usrmem )
130#ifdef dec_directives_win32
131!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: QP_ReadMatrix
132#endif
133 Use qp_data
134 implicit none
135 integer, intent (in) :: n ! number of variables
136 integer, intent (in) :: m ! number of constraints
137 integer, intent (in) :: nz ! number of nonzeros
138 real*8, intent (in out), dimension(n) :: lower ! vector of lower bounds
139 real*8, intent (in out), dimension(n) :: curr ! vector of initial values
140 real*8, intent (in out), dimension(n) :: upper ! vector of upper bounds
141 integer, intent (in out), dimension(n) :: vsta ! vector of initial variable status
142 ! (not defined here)
143 integer, intent (out), dimension(m) :: type ! vector of equation types
144 integer, intent (in out), dimension(m) :: esta ! vector of initial equation status
145 ! (not defined here)
146 real*8, intent (in out), dimension(m) :: rhs ! vector of right hand sides
147 integer, intent (in out), dimension(n+1) :: colsta ! vector with start of column indices
148 integer, intent (out), dimension(nz) :: rowno ! vector of row numbers
149 integer, intent (in out), dimension(nz) :: nlflag ! vector of nonlinearity flags
150 real*8, intent (in out), dimension(nz) :: value ! vector of matrix values
151 real*8 usrmem(*) ! optional user memory
152 integer :: i, j
153!
154! Information about Variables:
155! Default: Lower = -Inf, Curr = 0, and Upper = +inf.
156! Default: the status information in Vsta is not used.
157!
158! Nondefault: Lower bound = 0 for all variables
159!
160 do i = 1, nn+1
161 lower(i) = 0.0d0
162 enddo
163!
164! Information about Constraints:
165! Default: Rhs = 0
166! Default: the status information in Esta and the function
167! value in FV are not used.
168! Default: Type: There is no default.
169! 0 = Equality,
170! 1 = Greater than or equal,
171! 2 = Less than or equal,
172! 3 = Non binding.
173!
174! Constraint 1, defines Objective variable
175! Rhs = 0.0 and type Equality
176!
177 type(1) = 0
178!
179! Constraint 2 (Sum of all vars = 1)
180! Rhs = 1 and type Equality
181!
182 rhs(2) = 1.d0
183 type(2) = 0
184!
185! Information about the Jacobian. CONOPT expects a columnwise
186! representation in Rowno, Value, Nlflag and Colsta.
187!
188! Colsta = Start of column indices (No Defaults):
189! Rowno = Row indices
190! Value = Value of derivative (by default only linear
191! derivatives are used)
192! Nlflag = 0 for linear and 1 for nonlinear derivative
193! (not needed for completely linear models)
194!
195 j = 1 ! counter for current nonzero
196 do i = 1, nn
197 colsta(i) = j
198 rowno(j) = 1
199 nlflag(j) = 1
200 j = j + 1
201 rowno(j) = 2
202 value(j) = 1.d0
203 nlflag(j) = 0
204 j = j + 1
205 enddo
206 colsta(nn+1) = j
207 rowno(j) = 1
208 value(j) = -1.d0
209 nlflag(j) = 0
210 j = j + 1
211 colsta(nn+2) = j
212
213 qp_readmatrix = 0 ! Return value means OK
214
215end Function qp_readmatrix
216
217!> Compute nonlinear terms and non-constant Jacobian elements
218!!
219!! @include{doc} fdeval_params.dox
220Integer Function qp_fdeval( x, g, jac, rowno, jcnm, mode, ignerr, errcnt, &
221 n, nz, thread, usrmem )
222#ifdef dec_directives_win32
223!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: QP_FDEval
224#endif
225 Use qp_data
226 implicit none
227 integer, intent (in) :: n ! number of variables
228 integer, intent (in) :: rowno ! number of the row to be evaluated
229 integer, intent (in) :: nz ! number of nonzeros in this row
230 real*8, intent (in), dimension(n) :: x ! vector of current solution values
231 real*8, intent (in out) :: g ! constraint value
232 real*8, intent (in out), dimension(n) :: jac ! vector of derivatives for current constraint
233 integer, intent (in), dimension(nz) :: jcnm ! list of variables that appear nonlinearly
234 ! in this row. Ffor information only.
235 integer, intent (in) :: mode ! evaluation mode: 1 = function value
236 ! 2 = derivatives, 3 = both
237 integer, intent (in) :: ignerr ! if 1 then errors can be ignored as long
238 ! as errcnt is incremented
239 integer, intent (in out) :: errcnt ! error counter to be incremented in case
240 ! of function evaluation errors.
241 integer, intent (in) :: thread
242 real*8 usrmem(*) ! optional user memory
243 Integer :: i, j, k
244 real*8 :: sum
245!
246! Row 1: The objective function is nonlinear
247!
248 if ( rowno .eq. 1 ) then
249!
250! Mode = 1 or 3. Function value: G = x * Q * x / 2
251!
252 if ( mode .eq. 1 .or. mode .eq. 3 ) then
253 sum = 0.d0
254 do k = 1, nq
255 i = qrow(k)
256 j = qcol(k)
257 if ( i == j ) then
258 sum = sum + (x(i)-target(i))*q(k)*(x(j)-target(j))
259 else
260 sum = sum + 2*(x(i)-target(i))*q(k)*(x(j)-target(j))
261 endif
262 enddo
263 g = sum / 2.d0
264 endif
265!
266! Mode = 2 or 3: Derivative values: Q*x
267!
268 if ( mode .eq. 2 .or. mode .eq. 3 ) then
269 do i = 1, nn
270 jac(i) = 0
271 enddo
272 do k = 1, nq
273 i = qrow(k)
274 j = qcol(k)
275 if ( i == j ) then
276 jac(i) = jac(i) + q(k) * (x(i)-target(i))
277 else
278 jac(i) = jac(i) + q(k) * (x(j)-target(j))
279 jac(j) = jac(j) + q(k) * (x(i)-target(i))
280 endif
281 enddo
282 endif
283!
284! Row = 2: The row is linear and will not be called.
285!
286 endif
287 qp_fdeval = 0
288
289end Function qp_fdeval
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_optfile(cntvect, optfile)
define callback routine for defining an options file.
Definition conopt.f90:928
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_objvar(cntvect, objvar)
defines the Objective Variable.
Definition conopt.f90:257
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
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