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