CONOPT
Loading...
Searching...
No Matches
circles.f90
Go to the documentation of this file.
1!> @file circles.f90
2!! @ingroup FORT1THREAD_EXAMPLES
3!!
4!!
5!! Intersection of two circles with radius 1 and center in
6!! (0,0) and (1.6,0).
7!!
8!!
9!! For more information about the individual callbacks, please have a look at the source code.
10
11#if defined(_WIN32) && !defined(_WIN64)
12#define dec_directives_win32
13#endif
14
15!> Main program. A simple setup and call of CONOPT
16!!
17Program circles
18
20 Use conopt
21 IMPLICIT NONE
22!
23! Declare the user callback routines as Integer, External:
24!
25 Integer, External :: cir_readmatrix ! Mandatory Matrix definition routine defined below
26 Integer, External :: cir_fdeval ! Function and Derivative evaluation routine
27 ! needed a nonlinear model.
28 Integer, External :: std_status ! Standard callback for displaying solution status
29 Integer, External :: std_solution ! Standard callback for displaying solution values
30 Integer, External :: std_message ! Standard callback for managing messages
31 Integer, External :: std_errmsg ! Standard callback for managing error messages
32#ifdef dec_directives_win32
33!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Cir_ReadMatrix
34!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Cir_FDEval
35!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Status
36!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Solution
37!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Message
38!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_ErrMsg
39#endif
40!
41! Control vector
42!
43 INTEGER, Dimension(:), Pointer :: cntvect
44 INTEGER :: coi_error
45!
46! The equations are
47!
48! x**2 + y**2 = 1
49! (x-1.6)**2 + y**2 = 1
50!
51! starting in (0,1)
52!
53! Create and initialize a Control Vector
54!
55 call startup
56
57 coi_error = coi_create( cntvect )
58!
59! Tell CONOPT about the size of the model by populating the Control Vector:
60!
61 coi_error = max( coi_error, coidef_numvar( cntvect, 2 ) ) ! # variables
62 coi_error = max( coi_error, coidef_numcon( cntvect, 2 ) ) ! # constraints
63 coi_error = max( coi_error, coidef_numnz( cntvect, 4 ) ) ! # nonzeros in the Jacobian
64 coi_error = max( coi_error, coidef_numnlnz( cntvect, 4 ) ) ! # of which are nonlinear
65 coi_error = max( coi_error, coidef_optdir( cntvect, 0 ) ) ! No objective function but also not square system
66 coi_error = max( coi_error, coidef_optfile( cntvect, 'circles.opt' ) )
67 coi_error = max( coi_error, coidef_fvinclin( cntvect, 1 ) ) ! Function values inlcude linear terms
68!
69! Tell CONOPT about the callback routines:
70!
71 coi_error = max( coi_error, coidef_readmatrix( cntvect, cir_readmatrix ) )
72 coi_error = max( coi_error, coidef_fdeval( cntvect, cir_fdeval ) )
73 coi_error = max( coi_error, coidef_status( cntvect, std_status ) )
74 coi_error = max( coi_error, coidef_solution( cntvect, std_solution ) )
75 coi_error = max( coi_error, coidef_message( cntvect, std_message ) )
76 coi_error = max( coi_error, coidef_errmsg( cntvect, std_errmsg ) )
77
78#if defined(CONOPT_LICENSE_INT_1) && defined(CONOPT_LICENSE_INT_2) && defined(CONOPT_LICENSE_INT_3) && defined(CONOPT_LICENSE_TEXT)
79 coi_error = max( coi_error, coidef_license( cntvect, conopt_license_int_1, conopt_license_int_2, conopt_license_int_3, conopt_license_text) )
80#endif
81
82 If ( coi_error .ne. 0 ) THEN
83 write(*,*)
84 write(*,*) '**** Fatal Error while loading CONOPT Callback routines.'
85 write(*,*)
86 call flog( "Skipping Solve due to setup errors", 1 )
87 ENDIF
88!
89! Start CONOPT:
90!
91 coi_error = coi_solve( cntvect )
92 If ( coi_error /= 0 ) then
93 call flog( "Solve 1: Errors encountered during solution", 1 )
94 elseif ( stacalls == 0 .or. solcalls == 0 ) then
95 call flog( "Solve 1: Status or Solution routine was not called", 1 )
96 elseif ( sstat /= 1 .or. mstat /= 2 ) then
97 call flog( "Solve 1: Solver and Model Status was not as expected (1,2)", 1 )
98 endif
99
100 write(*,*)
101 write(*,*) 'End of Circles Equation example. Return code=',coi_error
102
103 if ( coi_free(cntvect) /= 0 ) call flog( "Error while freeing control vector",1)
104
105 call flog( "Successful Solve", 0 )
106
107End Program circles
108!
109! ============================================================================
110! Define information about the model:
111!
112
113!> Define information about the model
114!!
115!! @include{doc} readMatrix_params.dox
116Integer Function cir_readmatrix( lower, curr, upper, vsta, type, rhs, esta, &
117 colsta, rowno, value, nlflag, n, m, nz, &
118 usrmem )
119#ifdef dec_directives_win32
120!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Cir_ReadMatrix
121#endif
122 IMPLICIT NONE
123 integer, intent (in) :: n ! number of variables
124 integer, intent (in) :: m ! number of constraints
125 integer, intent (in) :: nz ! number of nonzeros
126 real*8, intent (in out), dimension(n) :: lower ! vector of lower bounds
127 real*8, intent (in out), dimension(n) :: curr ! vector of initial values
128 real*8, intent (in out), dimension(n) :: upper ! vector of upper bounds
129 integer, intent (in out), dimension(n) :: vsta ! vector of initial variable status
130 ! (not defined here)
131 integer, intent (out), dimension(m) :: type ! vector of equation types
132 integer, intent (in out), dimension(m) :: esta ! vector of initial equation status
133 ! (not defined here)
134 real*8, intent (in out), dimension(m) :: rhs ! vector of right hand sides
135 integer, intent (in out), dimension(n+1) :: colsta ! vector with start of column indices
136 integer, intent (out), dimension(nz) :: rowno ! vector of row numbers
137 integer, intent (in out), dimension(nz) :: nlflag ! vector of nonlinearity flags
138 real*8, intent (in out), dimension(nz) :: value ! vector of matrix values
139 real*8 usrmem(*) ! optional user memory
140
141!
142! Information about Variables:
143! Default: Lower = -Inf, Curr = 0, and Upper = +inf.
144! Default: the status information in Vsta is not used.
145!
146 curr(2) = 1.d0
147!
148! Information about Constraints:
149! Default: Rhs = 0
150!
151 rhs(1) = 1.0d0
152 rhs(2) = 1.0d0
153!
154! Default: the status information in Esta and the function
155! value in FV are not used.
156! Default: Type: There is no default.
157! 0 = Equality,
158! 1 = Greater than or equal,
159! 2 = Less than or equal,
160! 3 = Non binding.
161!
162 type(1) = 0
163 type(2) = 0
164!
165! Information about the Jacobian. We have to define Rowno, Value,
166! Nlflag and Colsta.
167!
168! Colsta = Start of column indices (No Defaults):
169! Rowno = Row indices
170! Value = Value of derivative (by default only linear
171! derivatives are used)
172! Nlflag = 0 for linear and 1 for nonlinear derivative
173! (not needed for completely linear models)
174!
175! Indices
176! x(1) x(2)
177! 1: 1 3
178! 2: 2 4
179!
180! Nonlinearity Structure: L = 0 are linear and NL = 1 are nonlinear
181! x(1) x(2)
182! 1: NL NL
183! 2: NL NL
184!
185! Value (Linear only) -- None
186!
187 colsta(1) = 1
188 colsta(2) = 3
189 colsta(3) = 5
190 rowno(1) = 1
191 rowno(2) = 2
192 rowno(3) = 1
193 rowno(4) = 2
194 nlflag(1) = 1
195 nlflag(2) = 1
196 nlflag(3) = 1
197 nlflag(4) = 1
199 cir_readmatrix = 0 ! Return value means OK
200
201end Function cir_readmatrix
202!
203!==========================================================================
204! Compute nonlinear terms and non-constant Jacobian elements
205!
206
207!> Compute nonlinear terms and non-constant Jacobian elements
208!!
209!! @include{doc} fdeval_params.dox
210Integer Function cir_fdeval( x, g, jac, rowno, jcnm, mode, ignerr, errcnt, &
211 n, nz, thread, usrmem )
212#ifdef dec_directives_win32
213!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Cir_FDEval
214#endif
215 IMPLICIT NONE
216 integer, intent (in) :: n ! number of variables
217 integer, intent (in) :: rowno ! number of the row to be evaluated
218 integer, intent (in) :: nz ! number of nonzeros in this row
219 real*8, intent (in), dimension(n) :: x ! vector of current solution values
220 real*8, intent (in out) :: g ! constraint value
221 real*8, intent (in out), dimension(n) :: jac ! vector of derivatives for current constraint
222 integer, intent (in), dimension(nz) :: jcnm ! list of variables that appear nonlinearly
223 ! in this row. Ffor information only.
224 integer, intent (in) :: mode ! evaluation mode: 1 = function value
225 ! 2 = derivatives, 3 = both
226 integer, intent (in) :: ignerr ! if 1 then errors can be ignored as long
227 ! as errcnt is incremented
228 integer, intent (in out) :: errcnt ! error counter to be incremented in case
229 ! of function evaluation errors.
230 integer, intent (in) :: thread
231 real*8 usrmem(*) ! optional user memory
232
233 cir_fdeval = 0
234!
235!
236 if ( rowno .eq. 1 ) then
237!
238! Row 1: x**2 + y**2 = 1
239! Mode = 1 or 3. Function value:
240!
241 if ( mode .eq. 1 .or. mode .eq. 3 ) then
242 g = x(1)**2 + x(2)**2
243 write(10,"('Row 1: X=',1p,2e20.10,' G =',1p,e20.10)") x, g
244 endif
245!
246! Mode = 2 or 3: Derivative values:
247!
248 if ( mode .eq. 2 .or. mode .eq. 3 ) then
249 jac(1) = 2.0d0*x(1)
250 jac(2) = 2.0d0*x(2)
251 write(10,"('Row 1: X=',1p,2e20.10,' Jac=',1p,2e20.10)") x, jac
252 endif
253 elseif ( rowno .eq. 2 ) then
254!
255! Row 2: (x-1.6)**2 + y**2 = 1
256! Mode = 1 or 3. Function value:
257!
258 if ( mode .eq. 1 .or. mode .eq. 3 ) then
259 g = (x(1)-1.6d0)**2 + x(2)**2
260 write(10,"('Row 2: X=',1p,2e20.10,' G =',1p,e20.10)") x, g
261 endif
262!
263! Mode = 2 or 3: Derivative values:
264!
265 if ( mode .eq. 2 .or. mode .eq. 3 ) then
266 jac(1) = 2.0d0*(x(1)-1.6d0)
267 jac(2) = 2.0d0*x(2)
268 write(10,"('Row 2: X=',1p,2e20.10,' Jac=',1p,2e20.10)") x, jac
269 endif
270 else
271 cir_fdeval = 1
272 endif
273
274end Function cir_fdeval
program circles
Main program. A simple setup and call of CONOPT.
Definition circles.f90:19
integer function cir_fdeval(x, g, jac, rowno, jcnm, mode, ignerr, errcnt, n, nz, thread, usrmem)
Compute nonlinear terms and non-constant Jacobian elements.
Definition circles.f90:200
integer function cir_readmatrix(lower, curr, upper, vsta, type, rhs, esta, colsta, rowno, value, nlflag, n, m, nz, usrmem)
Define information about the model.
Definition circles.f90:110
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
integer function std_message(smsg, dmsg, nmsg, llen, usrmem, msgv)
Definition comdecl.f90:205
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_license(cntvect, licint1, licint2, licint3, licstring)
define the License Information.
Definition conopt.f90:293
integer(c_int) function coidef_fvinclin(cntvect, fvinclin)
include the linear terms in function evaluations.
Definition conopt.f90:1053
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 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
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
integer mstat
Definition comdecl.f90:17
subroutine startup
Definition comdecl.f90:41