CONOPT
Loading...
Searching...
No Matches
lp10.f90
Go to the documentation of this file.
1!> @file lp10.f90
2!! @ingroup FORT1THREAD_EXAMPLES
3!!
4!! Test model lp10.gms -- case with an error
5!!
6!! For more information about the individual callbacks, please have a look at the source code.
7
8
10 real*8 :: case_c
11end module caseconst
12!> Main program. A simple setup and call of CONOPT
13!!
14Program lp10
15
16 Use proginfo
17 Use coidef
18 Use caseconst
19 implicit None
20!
21! Declare the user callback routines as Integer, External:
22!
23 Integer, External :: lp10_readmatrix ! Mandatory Matrix definition routine defined below
24 Integer, External :: std_status ! Standard callback for displaying solution status
25 Integer, External :: std_solution ! Standard callback for displaying solution values
26 Integer, External :: std_message ! Standard callback for managing messages
27 Integer, External :: std_errmsg ! Standard callback for managing error messages
28#if defined(itl)
29!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lp10_ReadMatrix
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 Integer :: iter
43!
44! Create and initialize a Control Vector
45!
46 call startup
47
48 numcallback = coidef_size()
49 Allocate( cntvect(numcallback) )
50 coi_error = coidef_inifort( cntvect )
51!
52! Tell CONOPT about the size of the model by populating the Control Vector:
53!
54 coi_error = max( coi_error, coidef_numvar( cntvect, 5 ) ) ! # variables
55 coi_error = max( coi_error, coidef_numcon( cntvect, 2 ) ) ! # constraints
56 coi_error = max( coi_error, coidef_numnz( cntvect, 7 ) ) ! # nonzeros in the Jacobian
57 coi_error = max( coi_error, coidef_numnlnz( cntvect, 0 ) ) ! # of which are nonlinear
58 coi_error = max( coi_error, coidef_optdir( cntvect, -1 ) ) ! Minimize
59 coi_error = max( coi_error, coidef_objvar( cntvect, 5 ) ) ! Objective is variable 5
60 coi_error = max( coi_error, coidef_optfile( cntvect, 'lp10.opt' ) )
61!
62! Tell CONOPT about the callback routines:
63!
64 coi_error = max( coi_error, coidef_readmatrix( cntvect, lp10_readmatrix ) )
65 coi_error = max( coi_error, coidef_status( cntvect, std_status ) )
66 coi_error = max( coi_error, coidef_solution( cntvect, std_solution ) )
67 coi_error = max( coi_error, coidef_message( cntvect, std_message ) )
68 coi_error = max( coi_error, coidef_errmsg( cntvect, std_errmsg ) )
69
70#if defined(LICENSE_INT_1) && defined(LICENSE_INT_2) && defined(LICENSE_INT_3) && defined(LICENSE_TEXT)
71 coi_error = max( coi_error, coidef_license( cntvect, license_int_1, license_int_2, license_int_3, license_text) )
72#endif
73
74 If ( coi_error .ne. 0 ) THEN
75 write(*,*)
76 write(*,*) '**** Fatal Error while loading CONOPT4 Callback routines.'
77 write(*,*)
78 call flog( "Skipping Solve due to setup errors", 1 )
79 ENDIF
80!
81! Save the solution so we can check the duals:
82!
83 do_allocate = .true.
84!
85! Start CONOPT:
86!
87 case_c = 1.0d0
88 do iter = 1, 100
89 case_c = case_c / 10.d0
90 coi_error = coi_solve( cntvect )
91 If ( coi_error /= 0 ) then
92 call flog( "Errors encountered during solution", 1 )
93 elseif ( stacalls == 0 .or. solcalls == 0 ) then
94 call flog( "Status or Solution routine was not called", 1 )
95 elseif ( sstat /= 1 ) then
96 call flog( "Solver Status was 1 not as expected.", 1 )
97 else
98 if ( mstat == 3 ) then ! Unbounded
99 if ( c_unbnd == 0 ) then
100 call flog( "Unbounded model but count of unbounded variables was zero.", 1 )
101 endif
102 else if ( mstat == 1 ) then ! Optimal
103 if ( c_unbnd /= 0 ) then
104 call flog( "Optimal model but Count of unbounded variables was not zero.", 1 )
105 Endif
106 if ( case_c > 1.d-6 ) then
107 call flog( "Optimal model but case_c > 1.d-6 should give unbounded.", 1 )
108 endif
109 exit
110 else
111 call flog( "Model Status was not as expected 3 or 1.", 1 )
112 endif
113 endif
114 enddo
115
116 write(*,*)
117 write(*,*) 'End of lp10 example. Return code=',coi_error
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 lp10
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 lp10_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 :: Lp10_ReadMatrix
137#endif
138 use caseconst
139 implicit none
140 integer, intent (in) :: n ! number of variables
141 integer, intent (in) :: m ! number of constraints
142 integer, intent (in) :: nz ! number of nonzeros
143 real*8, intent (in out), dimension(n) :: lower ! vector of lower bounds
144 real*8, intent (in out), dimension(n) :: curr ! vector of initial values
145 real*8, intent (in out), dimension(n) :: upper ! vector of upper bounds
146 integer, intent (in out), dimension(n) :: vsta ! vector of initial variable status
147 ! (not defined here)
148 integer, intent (out), dimension(m) :: type ! vector of equation types
149 integer, intent (in out), dimension(m) :: esta ! vector of initial equation status
150 ! (not defined here)
151 real*8, intent (in out), dimension(m) :: rhs ! vector of right hand sides
152 integer, intent (in out), dimension(n+1) :: colsta ! vector with start of column indices
153 integer, intent (out), dimension(nz) :: rowno ! vector of row numbers
154 integer, intent (in out), dimension(nz) :: nlflag ! vector of nonlinearity flags
155 real*8, intent (in out), dimension(nz) :: value ! vector of matrix values
156 real*8 usrmem(*) ! optional user memory
157!
158! Information about Variables:
159! Default: Lower = -Inf, Curr = 0, and Upper = +inf.
160! Default: the status information in Vsta is not used.
161!
162 lower(3) = 0.0d0
163 lower(4) = 0.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) = 0
176 type(2) = 0
177!
178! Information about the Jacobian. We use the standard method with
179! Rowno, Value, Nlflag and Colsta and we do not use Colno.
180!
181! Colsta = Start of column indices (No Defaults):
182! Rowno = Row indices
183! Value = Value of derivative (by default only linear
184! derivatives are used)
185! Nlflag = 0 for linear and 1 for nonlinear derivative
186! (not needed for completely linear models)
187!
188! Indices
189! x1 x2 x3 x4 x5
190! 1: 1 4 7
191! 2: 2 3 5 6
192!
193 colsta(1) = 1
194 colsta(2) = 3
195 colsta(3) = 4
196 colsta(4) = 6
197 colsta(5) = 7
198 colsta(6) = 8
199 rowno(1) = 1
200 rowno(2) = 2
201 rowno(3) = 2
202 rowno(4) = 1
203 rowno(5) = 2
204 rowno(6) = 2
205 rowno(7) = 1
206!
207! Nonlinearity Structure: Not needed for a linear model
208!
209! Value (Linear only)
210! x1 x2 x3 x4 x5
211! 1: -1 -1 1
212! 2: 1 c 1 -1
213!
214 value(1) = -1.d0
215 value(2) = 1.d0
216 value(3) = case_c
217 value(4) = -1.d0
218 value(5) = 1.d0
219 value(6) = -1.d0
220 value(7) = 1.d0
221
222 lp10_readmatrix = 0 ! Return value means OK
223
224end Function lp10_readmatrix
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_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 function lp10_readmatrix(lower, curr, upper, vsta, type, rhs, esta, colsta, rowno, value, nlflag, n, m, nz, usrmem)
Define information about the model.
Definition lp10.f90:135
program lp10
Main program. A simple setup and call of CONOPT.
Definition lp10.f90:14
real *8 case_c
Definition lp10.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 mstat
Definition comdecl.f90:11
subroutine startup
Definition comdecl.f90:35
integer c_unbnd
Definition comdecl.f90:16