CONOPT
Loading...
Searching...
No Matches
circle.f90
Go to the documentation of this file.
1!> @file circle.f90
2!! @ingroup FORT1THREAD_EXAMPLES
3!!
4!!
5!! The model approximates the inner of a circle with center in `(0,0)` and
6!! radius 1 by NP tangents defines in NP points spread evenly on the
7!! circle
8!! In addition we add the constraint `y = 0.5` and maximize `x`.
9!!
10!!
11!! For more information about the individual callbacks, please have a look at the source code.
12
13Module points
14 Integer, Parameter :: np = 24
15 real*8, dimension(NP) :: xp, yp
16 real*8, Parameter :: pi = 3.141592d0
17end Module points
18
19!> Main program. A simple setup and call of CONOPT
20!!
21Program circle
22
23 Use proginfo
24 Use coidef
25 Use points
26 implicit None
27!
28! Declare the user callback routines as Integer, External:
29!
30 Integer, External :: circle_readmatrix ! Mandatory Matrix definition routine defined below
31 Integer, External :: circle_fdeval ! Function and Derivative evaluation routine
32 ! needed a nonlinear model.
33 Integer, External :: std_status ! Standard callback for displaying solution status
34 Integer, External :: std_solution ! Standard callback for displaying solution values
35 Integer, External :: std_message ! Standard callback for managing messages
36 Integer, External :: std_errmsg ! Standard callback for managing error messages
37 Integer, External :: std_triord ! Standard callback for Monongular order
38#if defined(itl)
39!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Circle_ReadMatrix
40!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Circle_FDEval
41!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Status
42!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Solution
43!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Message
44!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_ErrMsg
45!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_TriOrd
46#endif
47!
48! Control vector
49!
50 INTEGER :: numcallback
51 INTEGER, Dimension(:), Pointer :: cntvect
52 INTEGER :: coi_error
53 Integer :: i
54
55 call startup
56 do i = 1, np
57 xp(i) = sin(2*pi*i/np)
58 yp(i) = cos(2*pi*i/np)
59 enddo
60!
61! Create and initialize a Control Vector
62!
63 numcallback = coidef_size()
64 Allocate( cntvect(numcallback) )
65 coi_error = coidef_inifort( cntvect )
66!
67! Tell CONOPT about the size of the model by populating the Control Vector:
68!
69 coi_error = max( coi_error, coidef_numvar( cntvect, 2 ) ) ! # variables
70 coi_error = max( coi_error, coidef_numcon( cntvect, np ) ) ! # constraints
71 coi_error = max( coi_error, coidef_numnz( cntvect, 2*np ) ) ! # nonzeros in the Jacobian
72 coi_error = max( coi_error, coidef_numnlnz( cntvect, 0 ) ) ! # of which are nonlinear
73 coi_error = max( coi_error, coidef_optdir( cntvect, +1 ) ) ! Maximize
74 coi_error = max( coi_error, coidef_objvar( cntvect, 1 ) ) ! Objective is variable 3
75 coi_error = max( coi_error, coidef_optfile( cntvect, 'Circle.opt' ) )
76!
77! Tell CONOPT about the callback routines:
78!
79 coi_error = max( coi_error, coidef_readmatrix( cntvect, circle_readmatrix ) )
80 coi_error = max( coi_error, coidef_fdeval( cntvect, circle_fdeval ) )
81 coi_error = max( coi_error, coidef_status( cntvect, std_status ) )
82 coi_error = max( coi_error, coidef_solution( cntvect, std_solution ) )
83 coi_error = max( coi_error, coidef_message( cntvect, std_message ) )
84 coi_error = max( coi_error, coidef_errmsg( cntvect, std_errmsg ) )
85 coi_error = max( coi_error, coidef_triord( cntvect, std_triord ) )
86
87#if defined(LICENSE_INT_1) && defined(LICENSE_INT_2) && defined(LICENSE_INT_3) && defined(LICENSE_TEXT)
88 coi_error = max( coi_error, coidef_license( cntvect, license_int_1, license_int_2, license_int_3, license_text) )
89#endif
90
91 If ( coi_error .ne. 0 ) THEN
92 write(*,*)
93 write(*,*) '**** Fatal Error while loading CONOPT Callback routines.'
94 write(*,*)
95 call flog( "Skipping Solve due to setup errors", 1 )
96 ENDIF
97!
98! Save the solution so we can check the duals:
99!
100 do_allocate = .true.
101!
102! Start CONOPT:
103!
104 coi_error = coi_solve( cntvect )
105
106 write(*,*)
107 write(*,*) 'End of Circle example. 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 /= 1 ) then
114 call flog( "Solver and Model Status was not as expected (1,1)", 1 )
115! elseif ( abs( OBJ-exp(2.0d0) ) > 0.000001d0 ) then ! not yet known
116! call flog( "Incorrect objective returned", 1 )
117 Else
118 Call checkdual( 'Circle', maximize )
119 endif
120
121 if ( coi_free(cntvect) /= 0 ) call flog( "Error while freeing control vector",1)
122
123 call flog( "Successful Solve", 0 )
124
125End Program circle
126!
127! ============================================================================
128! Define information about the model:
129!
130
131!> Define information about the model
132!!
133!! @include{doc} readMatrix_params.dox
134Integer Function circle_readmatrix( lower, curr, upper, vsta, type, rhs, esta, &
135 colsta, rowno, value, nlflag, n, m, nz, &
136 usrmem )
137#if defined(itl)
138!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Circle_ReadMatrix
139#endif
140 Use points
141 implicit none
142 integer, intent (in) :: n ! number of variables
143 integer, intent (in) :: m ! number of constraints
144 integer, intent (in) :: nz ! number of nonzeros
145 real*8, intent (in out), dimension(n) :: lower ! vector of lower bounds
146 real*8, intent (in out), dimension(n) :: curr ! vector of initial values
147 real*8, intent (in out), dimension(n) :: upper ! vector of upper bounds
148 integer, intent (in out), dimension(n) :: vsta ! vector of initial variable status
149 ! (not defined here)
150 integer, intent (out), dimension(m) :: type ! vector of equation types
151 integer, intent (in out), dimension(m) :: esta ! vector of initial equation status
152 ! (not defined here)
153 real*8, intent (in out), dimension(m) :: rhs ! vector of right hand sides
154 integer, intent (in out), dimension(n+1) :: colsta ! vector with start of column indices
155 integer, intent (out), dimension(nz) :: rowno ! vector of row numbers
156 integer, intent (in out), dimension(nz) :: nlflag ! vector of nonlinearity flags
157 real*8, intent (in out), dimension(nz) :: value ! vector of matrix values
158 real*8 usrmem(*) ! optional user memory
159 Integer :: i
160!
161! Information about Variables:
162! Default: Lower = -Inf, Curr = 0, and Upper = +inf.
163! Default: the status information in Vsta is not used.
164!
165! The model uses defaults
166!
167! Information about Constraints:
168! Default: Rhs = 0
169! Default: the status information in Esta and the function
170! value in FV are not used.
171! Default: Type: There is no default.
172! 0 = Equality,
173! 1 = Greater than or equal,
174! 2 = Less than or equal,
175! 3 = Non binding.
176!
177! All Constraints:
178! Rhs = 1.0 and type Less than or Equal
179!
180 do i = 1, np
181 rhs(i) = 1.0d0
182 type(i) = 2
183 enddo
184 curr(1) = 1.0d0
185 lower(2) = 0.5d0
186 upper(2) = 0.5d0
187 curr(2) = 0.5d0
188!
189! Information about the Jacobian. We use the standard method with
190! Rowno, Value, Nlflag and Colsta and we do not use Colno.
191!
192! Colsta = Start of column indices (No Defaults):
193! Rowno = Row indices
194! Value = Value of derivative (by default only linear
195! derivatives are used)
196! Nlflag = 0 for linear and 1 for nonlinear derivative
197! (not needed for completely linear models)
198!
199! Indices
200! x(1) x(2)
201! 1: 1 NP+1
202! 2: 2 NP+2
203! ...
204! NP: NP 2*NP
205!
206 colsta(1) = 1
207 colsta(2) = np+1
208 colsta(3) = 2*np+1
209 do i = 1, np
210 rowno(i) = i
211 rowno(i+np) = i
212 enddo
213!
214! Nonlinearity Structure: L = 0 are linear and NL = 1 are nonlinear
215! The model is linear and nlflag is not needed
216!
217! Value (Linear only)
218! x(1) x(2)
219! i: xp(i) yp(i)
220!
221 do i = 1, np
222 value(i) = xp(i)
223 value(i+np) = yp(i)
224 enddo
225 circle_readmatrix = 0 ! Return value means OK
226
227end Function circle_readmatrix
228!
229!==========================================================================
230! Compute nonlinear terms and non-constant Jacobian elements
231!
232
233!> Compute nonlinear terms and non-constant Jacobian elements
234!!
235!! @include{doc} fdeval_params.dox
236Integer Function circle_fdeval( x, g, jac, rowno, jcnm, mode, ignerr, errcnt, &
237 n, nz, thread, usrmem )
238#if defined(itl)
239!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Circle_FDEval
240#endif
241 Use points
242 implicit none
243 integer, intent (in) :: n ! number of variables
244 integer, intent (in) :: rowno ! number of the row to be evaluated
245 integer, intent (in) :: nz ! number of nonzeros in this row
246 real*8, intent (in), dimension(n) :: x ! vector of current solution values
247 real*8, intent (in out) :: g ! constraint value
248 real*8, intent (in out), dimension(n) :: jac ! vector of derivatives for current constraint
249 integer, intent (in), dimension(nz) :: jcnm ! list of variables that appear nonlinearly
250 ! in this row. Ffor information only.
251 integer, intent (in) :: mode ! evaluation mode: 1 = function value
252 ! 2 = derivatives, 3 = both
253 integer, intent (in) :: ignerr ! if 1 then errors can be ignored as long
254 ! as errcnt is incremented
255 integer, intent (in out) :: errcnt ! error counter to be incremented in case
256 ! of function evaluation errors.
257 integer, intent (in) :: thread
258 real*8 usrmem(*) ! optional user memory
259!
260! There are no nonlinear rows:
261!
262 circle_fdeval = 1
263
264end Function circle_fdeval
265
integer function circle_fdeval(x, g, jac, rowno, jcnm, mode, ignerr, errcnt, n, nz, thread, usrmem)
Compute nonlinear terms and non-constant Jacobian elements.
Definition circle.f90:238
program circle
Main program. A simple setup and call of CONOPT.
Definition circle.f90:21
integer function circle_readmatrix(lower, curr, upper, vsta, type, rhs, esta, colsta, rowno, value, nlflag, n, m, nz, usrmem)
Define information about the model.
Definition circle.f90:137
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_triord(mode, type, status, irow, icol, inf, value, resid, usrmem)
Definition comdecl.f90:291
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_triord(cntvect, coi_triord)
define callback routine for providing the triangular order information.
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
real *8, dimension(np) xp
Definition circle.f90:15
integer, parameter np
Definition circle.f90:14
real *8, parameter pi
Definition circle.f90:16
real *8, dimension(np) yp
Definition circle.f90:15
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