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