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