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