CONOPT
Loading...
Searching...
No Matches
const08.f90
Go to the documentation of this file.
1!> @file const08.f90
2!! @ingroup FORT1THREAD_EXAMPLES
3!!
4!!
5!! Model with derivatives that become constant after other variables are fixed.
6!!
7!! This is a CONOPT implementation of the GAMS model:
8!!
9!! @verbatim
10!! e1: max x1+x3
11!! e2: x1*x2 + log(x3) =l= 2.5
12!! x2.fx = 1;
13!! 1 <= x1 <= 3; x1.l = 2
14!! 1 <= x3 <= 5; x3.l = 3
15!! @endverbatim
16!!
17!! e1 is post-triangular.
18!! e2 has a combination of a constant Jacobian element and a monotone term.
19!! The expected solution is x1 = 1 and x3 = exp(1.5)
20!! Using monotonicity the solution may be proven globally optimal.
21!!
22!!
23!! For more information about the individual callbacks, please have a look at the source code.
24
25!> Main program. A simple setup and call of CONOPT
26!!
27Program const08
28
29 Use proginfo
30 Use coidef
31 implicit None
32!
33! Declare the user callback routines as Integer, External:
34!
35 Integer, External :: con_readmatrix ! Mandatory Matrix definition routine defined below
36 Integer, External :: con_fdeval ! Function and Derivative evaluation routine
37 ! needed a nonlinear model.
38 Integer, External :: con_fdinterval ! Function and Derivative evaluation routine
39 ! optional for a nonlinear model.
40 Integer, External :: std_status ! Standard callback for displaying solution status
41 Integer, External :: std_solution ! Standard callback for displaying solution values
42 Integer, External :: std_message ! Standard callback for managing messages
43 Integer, External :: std_errmsg ! Standard callback for managing error messages
44#if defined(itl)
45!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Con_ReadMatrix
46!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Con_FDEval
47!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Con_FDInterval
48!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Status
49!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Solution
50!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Message
51!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_ErrMsg
52#endif
53!
54! Control vector
55!
56 INTEGER, Dimension(:), Pointer :: cntvect
57 INTEGER :: coi_error
58!
59! Create and initialize a Control Vector
60!
61 call startup
62
63 coi_error = coi_createfort( cntvect )
64!
65! Tell CONOPT about the size of the model by populating the Control Vector:
66!
67 coi_error = max( coi_error, coidef_numvar( cntvect, 3 ) ) ! # variables
68 coi_error = max( coi_error, coidef_numcon( cntvect, 2 ) ) ! # constraints
69 coi_error = max( coi_error, coidef_numnz( cntvect, 5 ) ) ! # nonzeros in the Jacobian
70 coi_error = max( coi_error, coidef_numnlnz( cntvect, 3 ) ) ! # of which are nonlinear
71 coi_error = max( coi_error, coidef_optdir( cntvect, 1 ) ) ! Maximize
72 coi_error = max( coi_error, coidef_objcon( cntvect, 1 ) ) ! Objective is constraint 1
73 coi_error = max( coi_error, coidef_optfile( cntvect, 'const08.opt' ) )
74!
75! Tell CONOPT about the callback routines:
76!
77 coi_error = max( coi_error, coidef_readmatrix( cntvect, con_readmatrix ) )
78 coi_error = max( coi_error, coidef_fdeval( cntvect, con_fdeval ) )
79 coi_error = max( coi_error, coidef_fdinterval( cntvect, con_fdinterval ) )
80 coi_error = max( coi_error, coidef_status( cntvect, std_status ) )
81 coi_error = max( coi_error, coidef_solution( cntvect, std_solution ) )
82 coi_error = max( coi_error, coidef_message( cntvect, std_message ) )
83 coi_error = max( coi_error, coidef_errmsg( cntvect, std_errmsg ) )
84
85#if defined(LICENSE_INT_1) && defined(LICENSE_INT_2) && defined(LICENSE_INT_3) && defined(LICENSE_TEXT)
86 coi_error = max( coi_error, coidef_license( cntvect, license_int_1, license_int_2, license_int_3, license_text) )
87#endif
88
89 If ( coi_error .ne. 0 ) THEN
90 write(*,*)
91 write(*,*) '**** Fatal Error while loading CONOPT Callback routines.'
92 write(*,*)
93 call flog( "Skipping Solve due to setup errors", 1 )
94 ENDIF
95!
96! Save the solution so we can check the duals:
97!
98 do_allocate = .true.
99!
100! Start CONOPT:
101!
102 coi_error = coi_solve( cntvect )
103
104 write(*,*)
105 write(*,*) 'End of const08 example. Return code=',coi_error
106
107 If ( coi_error /= 0 ) then
108 call flog( "Errors encountered during solution", 1 )
109 elseif ( stacalls == 0 .or. solcalls == 0 ) then
110 call flog( "Status or Solution routine was not called", 1 )
111 elseif ( sstat /= 1 .or. ( mstat /= 1 .and. mstat /= 2 ) ) then
112 call flog( "Solver and Model Status was not as expected (1,1) or (1,2)", 1 )
113 elseif ( abs( obj-(1.0d0+exp(1.5d0))) > 0.000001d0 ) then
114 call flog( "Incorrect objective returned", 1 )
115 Else
116 Call checkdual( 'const08', maximize )
117 endif
118
119 if ( coi_free(cntvect) /= 0 ) call flog( "Error while freeing control vector",1)
120
121 call flog( "Successful Solve", 0 )
122
123End Program const08
124!
125! ============================================================================
126! Define information about the model:
127!
128
129!> Define information about the model
130!!
131!! @include{doc} readMatrix_params.dox
132Integer Function con_readmatrix( lower, curr, upper, vsta, type, rhs, esta, &
133 colsta, rowno, value, nlflag, n, m, nz, &
134 usrmem )
135#if defined(itl)
136!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Con_ReadMatrix
137#endif
138 implicit none
139 integer, intent (in) :: n ! number of variables
140 integer, intent (in) :: m ! number of constraints
141 integer, intent (in) :: nz ! number of nonzeros
142 real*8, intent (in out), dimension(n) :: lower ! vector of lower bounds
143 real*8, intent (in out), dimension(n) :: curr ! vector of initial values
144 real*8, intent (in out), dimension(n) :: upper ! vector of upper bounds
145 integer, intent (in out), dimension(n) :: vsta ! vector of initial variable status
146 ! (not defined here)
147 integer, intent (out), dimension(m) :: type ! vector of equation types
148 integer, intent (in out), dimension(m) :: esta ! vector of initial equation status
149 ! (not defined here)
150 real*8, intent (in out), dimension(m) :: rhs ! vector of right hand sides
151 integer, intent (in out), dimension(n+1) :: colsta ! vector with start of column indices
152 integer, intent (out), dimension(nz) :: rowno ! vector of row numbers
153 integer, intent (in out), dimension(nz) :: nlflag ! vector of nonlinearity flags
154 real*8, intent (in out), dimension(nz) :: value ! vector of matrix values
155 real*8 usrmem(*) ! optional user memory
156!
157! Information about Variables:
158! Default: Lower = -Inf, Curr = 0, and Upper = +inf.
159! Default: the status information in Vsta is not used.
160!
161 lower(1) = 1.0d0; curr(1) = 2.0d0; upper(1) = 3.0d0
162 lower(2) = 1.0d0; curr(2) = 1.0d0; upper(2) = 1.0d0
163 lower(3) = 1.0d0; curr(3) = 3.0d0; upper(3) = 5.0d0
164!
165! Information about Constraints:
166! Default: Rhs = 0
167! Default: the status information in Esta and the function
168! value in FV are not used.
169! Default: Type: There is no default.
170! 0 = Equality,
171! 1 = Greater than or equal,
172! 2 = Less than or equal,
173! 3 = Non binding.
174!
175 type(1) = 3
176 type(2) = 0
177 rhs(2) = 2.5d0
178!
179! Information about the Jacobian. We use the standard method with
180! Rowno, Value, Nlflag and Colsta and we do not use Colno.
181!
182! Colsta = Start of column indices (No Defaults):
183! Rowno = Row indices
184! Value = Value of derivative (by default only linear
185! derivatives are used)
186! Nlflag = 0 for linear and 1 for nonlinear derivative
187! (not needed for completely linear models)
188!
189! Indices
190! x(1) x(2) x(3)
191! 1: 1 4
192! 2: 2 3 5
193!
194 colsta(1) = 1
195 colsta(2) = 3
196 colsta(3) = 4
197 colsta(4) = 6
198 rowno(1) = 1
199 rowno(2) = 2
200 rowno(3) = 2
201 rowno(4) = 1
202 rowno(5) = 2
203!
204! Nonlinearity Structure: L = 0 are linear and NL = 1 are nonlinear
205! x(1) x(2) x(3)
206! 1: L L
207! 2: NL NL NL
208!
209 nlflag(1) = 0
210 nlflag(2) = 1
211 nlflag(3) = 1
212 nlflag(4) = 0
213 nlflag(5) = 1
214!
215! Value (Linear only)
216! x(1) x(2) x(3)
217! 1: +1 +1
218! 2: NL NL NL
219!
220 value(1) = +1.d0
221 value(4) = +1.d0
222
223 con_readmatrix = 0 ! Return value means OK
224
225end Function con_readmatrix
226!
227!==========================================================================
228! Compute nonlinear terms and non-constant Jacobian elements
229!
230
231!> Compute nonlinear terms and non-constant Jacobian elements
232!!
233!! @include{doc} fdeval_params.dox
234Integer Function con_fdeval( x, g, jac, rowno, jcnm, mode, ignerr, errcnt, &
235 n, nz, thread, usrmem )
236#if defined(itl)
237!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Con_FDEval
238#endif
239 implicit none
240 integer, intent (in) :: n ! number of variables
241 integer, intent (in) :: rowno ! number of the row to be evaluated
242 integer, intent (in) :: nz ! number of nonzeros in this row
243 real*8, intent (in), dimension(n) :: x ! vector of current solution values
244 real*8, intent (in out) :: g ! constraint value
245 real*8, intent (in out), dimension(n) :: jac ! vector of derivatives for current constraint
246 integer, intent (in), dimension(nz) :: jcnm ! list of variables that appear nonlinearly
247 ! in this row. Ffor information only.
248 integer, intent (in) :: mode ! evaluation mode: 1 = function value
249 ! 2 = derivatives, 3 = both
250 integer, intent (in) :: ignerr ! if 1 then errors can be ignored as long
251 ! as errcnt is incremented
252 integer, intent (in out) :: errcnt ! error counter to be incremented in case
253 ! of function evaluation errors.
254 integer, intent (in) :: thread
255 real*8 usrmem(*) ! optional user memory
256!
257! Row 1: the objective function is nonlinear
258!
259 if ( rowno .eq. 2 ) then
260!
261! Mode = 1 or 3: Function value
262!
263 if ( mode .eq. 1 .or. mode .eq. 3 ) then
264 g = x(1)*x(2) + log(x(3))
265 endif
266!
267! Mode = 2 or 3: Derivatives
268!
269 if ( mode .eq. 2 .or. mode .eq. 3 ) then
270 jac(1) = x(2)
271 jac(2) = x(1)
272 jac(3) = 1.0d0/x(3)
273 endif
274 con_fdeval = 0
275 Else
276 con_fdeval = 1 ! Should not happen
277 endif
278
279end Function con_fdeval
280
281
282!> Evaluating nonlinear functions and derivatives on an interval. Used in preprocessing
283!!
284!! @include{doc} fdinterval_params.dox
285Integer Function con_fdinterval( XMIN, XMAX, GMIN, GMAX, &
286 JMIN, JMAX, ROWNO, JCNM, &
287 MODE, PINF, N, NJ, USRMEM )
288#if defined(itl)
289!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Con_FDInterval
290#endif
291 Implicit None
292 INTEGER, Intent(IN) :: rowno, mode, n, nj
293 INTEGER, Dimension(NJ), Intent(IN) :: jcnm
294 real*8, Dimension(N), Intent(IN) :: xmin, xmax
295 real*8, Intent(IN OUT) :: gmin, gmax
296 real*8, Dimension(N), Intent(IN OUT) :: jmin, jmax
297 real*8, Intent(IN) :: pinf
298 real*8, Intent(IN OUT) :: usrmem(*)
299
300!
301! Row 2: x1*x2+log(x3) ! with known positive values
302!
303 if ( rowno .eq. 2 ) then
304!
305! Mode = 1 or 3. Function value
306!
307 if ( mode .eq. 1 .or. mode .eq. 3 ) then
308 gmin = xmin(1)*xmin(2) + log(xmin(3))
309 gmax = xmax(1)*xmax(2) + log(xmax(3))
310 endif
311!
312! Mode = 2 or 3: Derivative values:
313!
314 if ( mode .eq. 2 .or. mode .eq. 3 ) then
315 jmin(1) = xmin(2)
316 jmin(2) = xmin(1)
317 jmin(3) = 1.0d0/xmax(3)
318 jmax(1) = xmax(2)
319 jmax(2) = xmax(1)
320 jmax(3) = 1.d0/xmin(3)
321 endif
323 else
324!
325! There are no other rows:
326!
328 endif
329
330end Function con_fdinterval
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 con_readmatrix(lower, curr, upper, vsta, type, rhs, esta, colsta, rowno, value, nlflag, n, m, nz, usrmem)
Define information about the model.
Definition const01.f90:133
integer function con_fdeval(x, g, jac, rowno, jcnm, mode, ignerr, errcnt, n, nz, thread, usrmem)
Compute nonlinear terms and non-constant Jacobian elements.
Definition const01.f90:238
integer function con_fdinterval(xmin, xmax, gmin, gmax, jmin, jmax, rowno, jcnm, mode, pinf, n, nj, usrmem)
Evaluating nonlinear functions and derivatives on an interval. Used in preprocessing.
Definition const01.f90:291
program const08
Main program. A simple setup and call of CONOPT.
Definition const08.f90:27
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_fdinterval(cntvect, coi_fdinterval)
define callback routine for performing function and derivative evaluations on intervals.
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 coi_solve(cntvect)
method for starting the solving process of CONOPT.
Definition coistart.f90:14
#define nj
Definition mp_trans.c:46
real *8 obj
Definition comdecl.f90:10
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
logical do_allocate
Definition comdecl.f90:21
integer, parameter maximize
Definition comdecl.f90:25
integer mstat
Definition comdecl.f90:11
subroutine startup
Definition comdecl.f90:35