CONOPT
Loading...
Searching...
No Matches
polygon.f90
Go to the documentation of this file.
1!> @file polygon.f90
2!! @ingroup FORT1THREAD_EXAMPLES
3!!
4!!
5!! Polygon model from COPS test set.
6!!
7!!
8!! This is a CONOPT implementation of the GAMS model:
9!!
10!! @verbatim
11!! set i sides / i1*i%nv% /;
12!!
13!! alias (i,j);
14!!
15!! scalar pi a famous constant;
16!!
17!!
18!! positive variables
19!! r(i) polar radius (distance to fixed vertex)
20!! theta(i) polar angle (measured from fixed direction)
21!! variable polygon_area
22!!
23!! equations
24!! obj
25!! distance(i,j)
26!! ordered(i) ;
27!!
28!! obj.. polygon_area =e= 0.5*sum(j(i+1), r(i+1)*r(i)*sin(theta(i+1)-theta(i)));
29!!
30!! ordered(i+1).. theta(i) =l= theta(i+1);
31!!
32!! distance(i,j)$(ord(j) > ord(i))..
33!! sqr(r(i)) + sqr(r(j)) - 2*r(i)*r(j)*cos(theta(j)-theta(i)) =l= 1;
34!!
35!!
36!! pi = 2*arctan(inf);
37!!
38!! r.up(i) = 1;
39!! theta.up(i) = pi;
40!!
41!! r.fx('i%nv%') = 0;
42!! theta.fx('i%nv%') = pi;
43!!
44!! r.l(i) = 4*ord(i)*(card(i)+ 1- ord(i))/sqr(card(i)+1);
45!! theta.l(i) = pi*ord(i)/card(i);
46!!
47!! model polygon /all/;
48!!
49!! $if set workspace polygon.workspace = %workspace%;
50!!
51!! solve polygon using nlp maximizing polygon_area;
52!! @endverbatim
53!!
54!!
55!! For more information about the individual callbacks, please have a look at the source code.
56
57Module pdata
58 Integer, Parameter :: np = 100 ! Number of polygon points
59 Integer, Parameter :: nv = 2*np ! Number of variables
60 Integer, Parameter :: no = np-1 ! Number of Ordered constraints
61 Integer, Parameter :: nd = np*(np-1)/2 ! Number of distance constraints
62 Integer, dimension(ND) :: rowmapi, rowmapj ! Map from distance index to i and j
63 Integer, dimension(NP,NP) :: rowmapk
64end module pdata
65
66!> Main program. A simple setup and call of CONOPT
67!!
68Program tutorial
69
70 Use proginfo
71 Use coidef
72 Use pdata
73 implicit None
74!
75! Declare the user callback routines as Integer, External:
76!
77 Integer, External :: poly_readmatrix ! Mandatory Matrix definition routine defined below
78 Integer, External :: poly_fdeval ! Function and Derivative evaluation routine
79 ! needed a nonlinear model.
80 Integer, External :: std_status ! Standard callback for displaying solution status
81 Integer, External :: std_solution ! Standard callback for displaying solution values
82 Integer, External :: std_message ! Standard callback for managing messages
83 Integer, External :: std_errmsg ! Standard callback for managing error messages
84#if defined(itl)
85!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Poly_ReadMatrix
86!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Poly_FDEval
87!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Status
88!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Solution
89!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Message
90!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_ErrMsg
91#endif
92!
93! Control vector
94!
95 INTEGER :: numcallback
96 INTEGER, Dimension(:), Pointer :: cntvect
97 INTEGER :: coi_error
98! Create and initialize a Control Vector
99!
100 call startup
101
102 numcallback = coidef_size()
103 Allocate( cntvect(numcallback) )
104 coi_error = coidef_inifort( cntvect )
105!
106! Tell CONOPT about the size of the model by populating the Control Vector:
107!
108 coi_error = max( coi_error, coidef_numvar( cntvect, nv ) ) ! # variables
109 coi_error = max( coi_error, coidef_numcon( cntvect, no + nd + 1 ) ) ! # constraints
110 coi_error = max( coi_error, coidef_numnz( cntvect, 2*no+4*nd+nv ) ) ! # nonzeros in the Jacobian
111 coi_error = max( coi_error, coidef_numnlnz( cntvect, 4*nd+nv ) ) ! # of which are nonlinear
112 coi_error = max( coi_error, coidef_optdir( cntvect, +1 ) ) ! Maximize
113 coi_error = max( coi_error, coidef_objcon( cntvect, no + nd + 1 ) ) ! Objective is last constraint
114 coi_error = max( coi_error, coidef_optfile( cntvect, 'polygon.opt' ) )
115!
116! Tell CONOPT about the callback routines:
117!
118 coi_error = max( coi_error, coidef_readmatrix( cntvect, poly_readmatrix ) )
119 coi_error = max( coi_error, coidef_fdeval( cntvect, poly_fdeval ) )
120 coi_error = max( coi_error, coidef_status( cntvect, std_status ) )
121 coi_error = max( coi_error, coidef_solution( cntvect, std_solution ) )
122 coi_error = max( coi_error, coidef_message( cntvect, std_message ) )
123 coi_error = max( coi_error, coidef_errmsg( cntvect, std_errmsg ) )
124
125#if defined(LICENSE_INT_1) && defined(LICENSE_INT_2) && defined(LICENSE_INT_3) && defined(LICENSE_TEXT)
126 coi_error = max( coi_error, coidef_license( cntvect, license_int_1, license_int_2, license_int_3, license_text) )
127#endif
128
129 If ( coi_error .ne. 0 ) THEN
130 write(*,*)
131 write(*,*) '**** Fatal Error while loading CONOPT Callback routines.'
132 write(*,*)
133 call flog( "Skipping Solve due to setup errors", 1 )
134 ENDIF
135!
136! Start CONOPT:
137!
138 coi_error = coi_solve( cntvect )
139 If ( coi_error /= 0 ) then
140 call flog( "Errors encountered during solution", 1 )
141 elseif ( stacalls == 0 .or. solcalls == 0 ) then
142 call flog( "Status or Solution routine was not called", 1 )
143 elseif ( sstat /= 1 .or. mstat /= 2 ) then
144 call flog( "Solver and Model Status was not as expected (1,2)", 1 )
145 endif
146
147 write(*,*)
148 write(*,*) 'End of Polygon example. Return code=',coi_error
149
150 if ( coi_free(cntvect) /= 0 ) call flog( "Error while freeing control vector",1)
151
152 call flog( "Successful Solve", 0 )
153
154End Program tutorial
155!
156! ============================================================================
157! Define information about the model:
158!
159
160!> Define information about the model
161!!
162!! @include{doc} readMatrix_params.dox
163Integer Function poly_readmatrix( lower, curr, upper, vsta, type, rhs, esta, &
164 colsta, rowno, value, nlflag, n, m, nz, &
165 usrmem )
166#if defined(itl)
167!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Poly_ReadMatrix
168#endif
169 Use pdata
170 implicit none
171 integer, intent (in) :: n ! number of variables
172 integer, intent (in) :: m ! number of constraints
173 integer, intent (in) :: nz ! number of nonzeros
174 real*8, intent (in out), dimension(n) :: lower ! vector of lower bounds
175 real*8, intent (in out), dimension(n) :: curr ! vector of initial values
176 real*8, intent (in out), dimension(n) :: upper ! vector of upper bounds
177 integer, intent (in out), dimension(n) :: vsta ! vector of initial variable status
178 ! (not defined here)
179 integer, intent (out), dimension(m) :: type ! vector of equation types
180 integer, intent (in out), dimension(m) :: esta ! vector of initial equation status
181 ! (not defined here)
182 real*8, intent (in out), dimension(m) :: rhs ! vector of right hand sides
183 integer, intent (in out), dimension(n+1) :: colsta ! vector with start of column indices
184 integer, intent (out), dimension(nz) :: rowno ! vector of row numbers
185 integer, intent (in out), dimension(nz) :: nlflag ! vector of nonlinearity flags
186 real*8, intent (in out), dimension(nz) :: value ! vector of matrix values
187 real*8 usrmem(*) ! optional user memory
188
189 Integer :: i, j, k
190 real*8, parameter :: pi = 3.141592
191
192!
193! Information about Variables:
194! Default: Lower = -Inf, Curr = 0, and Upper = +inf.
195! Default: the status information in Vsta is not used.
196!
197! r.up(i) = 1;
198! theta.up(i) = pi;
199!
200! r.fx('i%nv%') = 0;
201! theta.fx('i%nv%') = pi;
202!
203! r.l(i) = 4*ord(i)*(card(i)+ 1- ord(i))/sqr(card(i)+1);
204! theta.l(i) = pi*ord(i)/card(i);
205 DO i = 1, np-1
206 lower(i) = 0.d0
207 lower(i+np) = 0.d0
208 upper(i) = 1.d0
209 upper(i+np) = pi
210 curr(i) = 4.d0*i*(np+1-i)/(1+np)**2
211 curr(i+np) = pi*i/np
212 enddo
213 lower(np) = 0.d0
214 upper(np) = 0.d0
215 curr(np) = 0.d0
216 lower(2*np) = pi
217 upper(2*np) = pi
218 curr(2*np) = pi
219!
220! Information about Constraints:
221! Default: Rhs = 0
222! Default: the status information in Esta and the function
223! value in FV are not used.
224! Default: Type: There is no default.
225! 0 = Equality,
226! 1 = Greater than or equal,
227! 2 = Less than or equal,
228! 3 = Non binding.
229!
230! ordered(i+1).. theta(i) =l= theta(i+1);
231!
232! distance(i,j)$(ord(j) > ord(i))..
233! sqr(r(i)) + sqr(r(j)) - 2*r(i)*r(j)*cos(theta(j)-theta(i)) =l= 1;
234!
235! obj.. polygon_area =e= 0.5*sum(j(i+1), r(i+1)*r(i)*sin(theta(i+1)-theta(i)));
236!
237 k = 0
238 DO i = 1, no
239 k = k + 1
240 Type(k) = 2 ! Ordered constraint Less than and RHS = 0
241 Enddo
242 rowmapk = 0 ! Initialize so we can recognize (i,j) combinations that are not equations
243 DO i = 1, np-1
244 DO j = i+1, np
245 k = k + 1
246 Type(k) = 2 ! Dist(i,j) Less than
247 rhs(k) = 1.d0
248 rowmapi(k-no) = i ! To be used for easier lookup
249 rowmapj(k-no) = j
250 rowmapk(i,j) = k
251 Enddo
252 Enddo
253 k = k + 1
254 type(k) = 3 ! Objective of type Non binding
255!
256! Information about the Jacobian. We have to define Rowno, Value,
257! Nlflag and Colsta.
258!
259! Colsta = Start of column indices (No Defaults):
260! Rowno = Row indices
261! Value = Value of derivative (by default only linear
262! derivatives are used)
263! Nlflag = 0 for linear and 1 for nonlinear derivative
264! (not needed for completely linear models)
265!
266! Indices
267!
268! r(1) r(2) r(3) r(4) .. Th(1) th(2) th(3) th(4)
269! 1: 1 -1
270! 2: 1 -1
271! 3: 1 -1
272! ..
273! NP: NL NL NL NL
274! NP+1: NL NL NL NL
275! NP+2 NL NL NL NL
276! ..
277! NL NL NL NL
278! NL NL NL NL
279! ..
280! Obj: NL NL NL NL NL NL
281!
282! Nonlinearity Structure: L = 0 are linear and NL = 1 are nonlinear
283! All nonzeros are nonlinear
284!
285 k = 0
286 DO i = 1, np ! r variables
287 colsta(i) = k+1
288 DO j = i+1, np
289 k = k + 1
290 rowno(k) = rowmapk(i,j) ! Distance - i index
291 nlflag(k) = 1
292 enddo
293 DO j = 1,np
294 IF ( rowmapk(j,i) > 0 ) Then
295 k = k + 1
296 rowno(k) = rowmapk(j,i) ! Distance - j index
297 nlflag(k) = 1
298 endif
299 enddo
300 k = k + 1
301 rowno(k) = no+nd+1 ! Objective constraint
302 nlflag(k) = 1
303 enddo
304 DO i = 1, np ! theta variables
305 colsta(np+i) = k+1
306 if ( i < np ) then
307 k = k + 1
308 rowno(k) = i ! ordered constraint index i
309 nlflag(k) = 0
310 value(k) = 1.d0
311 endif
312 if ( i > 1 ) then
313 k = k + 1
314 rowno(k) = i-1 ! ordered constraint index i+1
315 nlflag(k) = 0
316 value(k) = -1.d0
317 endif
318 DO j = i+1, np
319 k = k + 1
320 rowno(k) = rowmapk(i,j) ! Distance - i index
321 nlflag(k) = 1
322 enddo
323 DO j = 1,np
324 IF ( rowmapk(j,i) > 0 ) Then
325 k = k + 1
326 rowno(k) = rowmapk(j,i) ! Distance - j index
327 nlflag(k) = 1
328 endif
329 enddo
330 k = k + 1
331 rowno(k) = no+nd+1 ! Objective constraint
332 nlflag(k) = 1
333 enddo
334 colsta(nv+1) = k+1
335
336 poly_readmatrix = 0 ! Return value means OK
337
338end Function poly_readmatrix
339!
340!==========================================================================
341! Compute nonlinear terms and non-constant Jacobian elements
342!
343
344!> Compute nonlinear terms and non-constant Jacobian elements
345!!
346!! @include{doc} fdeval_params.dox
347Integer Function poly_fdeval( x, g, jac, rowno, jcnm, mode, ignerr, errcnt, &
348 n, nz, thread, usrmem )
349#if defined(itl)
350!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Poly_FDEval
351#endif
352 Use pdata
353 implicit none
354 integer, intent (in) :: n ! number of variables
355 integer, intent (in) :: rowno ! number of the row to be evaluated
356 integer, intent (in) :: nz ! number of nonzeros in this row
357 real*8, intent (in), dimension(n) :: x ! vector of current solution values
358 real*8, intent (in out) :: g ! constraint value
359 real*8, intent (in out), dimension(n) :: jac ! vector of derivatives for current constraint
360 integer, intent (in), dimension(nz) :: jcnm ! list of variables that appear nonlinearly
361 ! in this row. Ffor information only.
362 integer, intent (in) :: mode ! evaluation mode: 1 = function value
363 ! 2 = derivatives, 3 = both
364 integer, intent (in) :: ignerr ! if 1 then errors can be ignored as long
365 ! as errcnt is incremented
366 integer, intent (in out) :: errcnt ! error counter to be incremented in case
367 ! of function evaluation errors.
368 integer, intent (in) :: thread
369 real*8 usrmem(*) ! optional user memory
370
371 Integer :: ne
372 Integer :: i, j
373 real*8 :: s, c
374
375 ne = n / 3
376!
377! Row NO+ND+1: the objective function is nonlinear
378!
379 if ( rowno == no+nd+1 ) then
380!
381! Mode = 1 or 3. Function value:
382! obj.. polygon_area =e= 0.5*sum(j(i+1), r(i+1)*r(i)*sin(theta(i+1)-theta(i)));
383!
384 if ( mode .eq. 1 .or. mode .eq. 3 ) then
385 g = 0.d0
386 do i = 1, np-1
387 g = g + x(i+1)*x(i)*sin(x(np+i+1)-x(np+i))
388 enddo
389 g = g * 0.5d0
390 endif
391!
392! Mode = 2 or 3: Derivative values: w.r.t. x: 2*(x(i)-x(j))*(-1/2)*dist(-3/2)
393!
394 if ( mode .eq. 2 .or. mode .eq. 3 ) then
395 do i = 1, nv
396 jac(i) = 0.d0
397 enddo
398 do i = 1, np-1
399! g = g + 0.5* x(i+1)*x(i)*Sin(x(NP+i+1)-x(NP+i))
400 s = 0.5d0*sin(x(np+i+1)-x(np+i))
401 jac(i+1) = jac(i+1) + x(i)*s
402 jac(i) = jac(i) + x(i+1)*s
403 c = 0.5d0*x(i+1)*x(i)*cos(x(np+i+1)-x(np+i))
404 jac(np+i+1) = jac(np+i+1) + c
405 jac(np+i) = jac(np+i) - c
406 enddo
407 endif
408!
409 Else if ( rowno > no ) then ! this is a distance constraint
410!
411! Mode = 1 or 3. Function value:
412! sqr(r(i)) + sqr(r(j)) - 2*r(i)*r(j)*cos(theta(j)-theta(i)) =l= 1;
413!
414 i = rowmapi(rowno-no)
415 j = rowmapj(rowno-no)
416 if ( mode .eq. 1 .or. mode .eq. 3 ) then
417 g = x(i)**2 + x(j)**2 -2.d0*x(i)*x(j)*cos(x(j+np)-x(i+np))
418 endif
419!
420! Mode = 2 or 3: Derivative values:
421!
422 if ( mode .eq. 2 .or. mode .eq. 3 ) then
423 c = cos(x(j+np)-x(i+np))
424 jac(i) = 2.d0*(x(i)-x(j)*c)
425 jac(j) = 2.d0*(x(j)-x(i)*c)
426 s = 2.d0*sin(x(j+np)-x(i+np))*x(i)*x(j)
427 jac(i+np) = -s
428 jac(j+np) = s
429 endif
430 Else
431 write(*,*) 'Fdeval called with rowno=',rowno
432 poly_fdeval = 1; return ! error in row number
433 endif
434 poly_fdeval = 0
435
436end Function poly_fdeval
Main program. A simple setup and call of CONOPT.
Definition tutorial.java:14
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_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_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 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, dimension(nd) rowmapj
Definition polygon.f90:62
integer, dimension(np, np) rowmapk
Definition polygon.f90:63
integer, dimension(nd) rowmapi
Definition polygon.f90:62
integer, parameter np
Definition polygon.f90:58
integer, parameter nd
Definition polygon.f90:61
integer, parameter nv
Definition polygon.f90:59
integer, parameter no
Definition polygon.f90:60
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
integer mstat
Definition comdecl.f90:11
subroutine startup
Definition comdecl.f90:35
integer function poly_readmatrix(lower, curr, upper, vsta, type, rhs, esta, colsta, rowno, value, nlflag, n, m, nz, usrmem)
Define information about the model.
Definition polygon.f90:166
integer function poly_fdeval(x, g, jac, rowno, jcnm, mode, ignerr, errcnt, n, nz, thread, usrmem)
Compute nonlinear terms and non-constant Jacobian elements.
Definition polygon.f90:349