CONOPT
Loading...
Searching...
No Matches
mp_polygon2.f90
Go to the documentation of this file.
1!> @file mp_polygon2.f90
2!! @ingroup FORTOPENMP_EXAMPLES
3!!
4!!
5!! Polygon model from COPS test set.
6!! The model is run sequentially and with various number of threads but with ThreadC defined.
7!! The models with various number of threads should produce identical solutions and we test
8!! this.
9!!
10!!
11!! For more information about the individual callbacks, please have a look at the source code.
12
13Module pdata2
14 Integer, Parameter :: np = 50 ! Number of polygon points
15 Integer, Parameter :: nv = 2*np ! Number of variables
16 Integer, Parameter :: no = np-1 ! Number of Ordered constraints
17 Integer, Parameter :: nd = np*(np-1)/2 ! Number of distance constraints
18 Integer, dimension(ND) :: rowmapi, rowmapj ! Map from distance index to i and j
19 Integer, dimension(NP,NP) :: rowmapk
20 real*8, dimension(NV) :: xprim, xprim1, xprim2, xprim4
21end module pdata2
22
23!> Main program. A simple setup and call of CONOPT
24!!
25Program polygon
26
27 Use proginfop
28 Use coidef
29 Use pdata2
30 Use omp_lib
31 implicit None
32!
33! Declare the user callback routines as Integer, External:
34!
35 Integer, External :: poly_readmatrix ! Mandatory Matrix definition routine defined below
36 Integer, External :: poly_fdeval ! Function and Derivative evaluation routine
37 ! needed a nonlinear model.
38 Integer, External :: std_status ! Standard callback for displaying solution status
39 Integer, External :: poly_solution ! Special callback for displaying and testing solution values
40 Integer, External :: std_message ! Standard callback for managing messages
41 Integer, External :: std_errmsg ! Standard callback for managing error messages
42#if defined(itl)
43!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Poly_ReadMatrix
44!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Poly_FDEval
45!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Status
46!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Poly_Solution
47!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Message
48!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_ErrMsg
49#endif
50!
51! Control vector
52!
53 INTEGER, Dimension(:), Pointer :: cntvect
54 INTEGER :: coi_error
55
56 real*8 time0, time1, time2, time4
57 Integer :: iter1, iter2, iter4
58 Integer :: i
59 Logical :: same12, same14, same24
60!
61! The model shown here in GAMS format
62!
63! set i sides / i1*i%nv% /;
64!
65! alias (i,j);
66!
67! scalar pi a famous constant;
68!
69!
70! positive variables
71! r(i) polar radius (distance to fixed vertex)
72! theta(i) polar angle (measured from fixed direction)
73! variable polygon_area
74!
75! equations
76! obj
77! distance(i,j)
78! ordered(i) ;
79!
80! obj.. polygon_area =e= 0.5*sum(j(i+1), r(i+1)*r(i)*sin(theta(i+1)-theta(i)));
81!
82! ordered(i+1).. theta(i) =l= theta(i+1);
83!
84! distance(i,j)$(ord(j) > ord(i))..
85! sqr(r(i)) + sqr(r(j)) - 2*r(i)*r(j)*cos(theta(j)-theta(i)) =l= 1;
86!
87!
88! pi = 2*arctan(inf);
89!
90! r.up(i) = 1;
91! theta.up(i) = pi;
92!
93! r.fx('i%nv%') = 0;
94! theta.fx('i%nv%') = pi;
95!
96! r.l(i) = 4*ord(i)*(card(i)+ 1- ord(i))/sqr(card(i)+1);
97! theta.l(i) = pi*ord(i)/card(i);
98!
99! model polygon /all/;
100!
101! $if set workspace polygon.workspace = %workspace%;
102!
103! solve polygon using nlp maximizing polygon_area;
104!
105! Create and initialize a Control Vector
106!
107 call startup
108
109 coi_error = coi_createfort( cntvect )
110!
111! Tell CONOPT about the size of the model by populating the Control Vector:
112!
113 coi_error = max( coi_error, coidef_numvar( cntvect, nv ) ) ! # variables
114 coi_error = max( coi_error, coidef_numcon( cntvect, no + nd + 1 ) ) ! # constraints
115 coi_error = max( coi_error, coidef_numnz( cntvect, 2*no+4*nd+nv ) ) ! # nonzeros in the Jacobian
116 coi_error = max( coi_error, coidef_numnlnz( cntvect, 4*nd+nv ) ) ! # of which are nonlinear
117 coi_error = max( coi_error, coidef_optdir( cntvect, +1 ) ) ! Maximize
118 coi_error = max( coi_error, coidef_objcon( cntvect, no + nd + 1 ) ) ! Objective is last constraint
119 coi_error = max( coi_error, coidef_optfile( cntvect, 'mp_polygon2.opt' ) )
120!
121! Tell CONOPT about the callback routines:
122!
123 coi_error = max( coi_error, coidef_readmatrix( cntvect, poly_readmatrix ) )
124 coi_error = max( coi_error, coidef_fdeval( cntvect, poly_fdeval ) )
125 coi_error = max( coi_error, coidef_status( cntvect, std_status ) )
126 coi_error = max( coi_error, coidef_solution( cntvect, poly_solution ) )
127 coi_error = max( coi_error, coidef_message( cntvect, std_message ) )
128 coi_error = max( coi_error, coidef_errmsg( cntvect, std_errmsg ) )
129
130#if defined(LICENSE_INT_1) && defined(LICENSE_INT_2) && defined(LICENSE_INT_3) && defined(LICENSE_TEXT)
131 coi_error = max( coi_error, coidef_license( cntvect, license_int_1, license_int_2, license_int_3, license_text) )
132#endif
133
134 If ( coi_error .ne. 0 ) THEN
135 write(*,*)
136 write(*,*) '**** Fatal Error while loading CONOPT Callback routines.'
137 write(*,*)
138 call flog( "Skipping Solve due to setup errors", 1 )
139 ENDIF
140!
141! Start CONOPT -- First time use ThreadC < ThreadS. We should get error 115.
142!
143 coi_error = max( coi_error, coidef_threads( cntvect, 2 ) ) !
144 coi_error = max( coi_error, coidef_threadc( cntvect, 1 ) ) !
145 coi_error = coi_solve( cntvect )
146
147 write(*,*)
148 write(*,*) 'End of Polygon example with one thread. Return code=',coi_error
149
150 If ( coi_error /= 115 ) then
151 call flog( "ThreadC = 1 and ThreadS = 2 did not give error return 115.", 1 )
152 endif
153 coi_error = max( coi_error, coidef_threadc( cntvect, 8 ) ) ! Cut work in 8 when testing threads
154!
155! Start CONOPT again -- this time exactly 1 thread
156!
157 coi_error = max( coi_error, coidef_threads( cntvect, 1 ) ) ! Use 1 thread
158 time0 = omp_get_wtime()
159 coi_error = coi_solve( cntvect )
160 time1 = omp_get_wtime() - time0
161
162 write(*,*)
163 write(*,*) 'End of Polygon example with 1 threads. Return code=',coi_error
164
165 If ( coi_error /= 0 ) then
166 call flog( "Single thread: Errors encountered during solution", 1 )
167 elseif ( stacalls == 0 .or. solcalls == 0 ) then
168 call flog( "Single thread: Status or Solution routine was not called", 1 )
169 elseif ( sstat /= 1 .or. mstat /= 2 ) then
170 call flog( "Single thread: Solver and Model Status was not as expected (1,2)", 1 )
171 endif
172 iter1 = miter
173 xprim1 = xprim
174!
175! Start CONOPT again -- this time exactly 2 threads -- we keep ThreadC = 8
176!
177 coi_error = max( coi_error, coidef_threads( cntvect, 2 ) ) ! Use 2 thread
178 time0 = omp_get_wtime()
179 coi_error = coi_solve( cntvect )
180 time2 = omp_get_wtime() - time0
181
182 write(*,*)
183 write(*,*) 'End of Polygon example with 2 threads. Return code=',coi_error
184
185 If ( coi_error /= 0 ) then
186 call flog( "Double threads: Errors encountered during solution", 1 )
187 elseif ( stacalls == 0 .or. solcalls == 0 ) then
188 call flog( "Double threads: Status or Solution routine was not called", 1 )
189 elseif ( sstat /= 1 .or. mstat /= 2 ) then
190 call flog( "Double threads: Solver and Model Status was not as expected (1,2)", 1 )
191 endif
192 iter2 = miter
193 xprim2 = xprim
194!
195! Start CONOPT again -- this time exactly 4 threads -- we keep ThreadC = 16
196!
197 coi_error = max( coi_error, coidef_threads( cntvect, 4 ) ) ! Use 4 thread
198 time0 = omp_get_wtime()
199 coi_error = coi_solve( cntvect )
200 time4 = omp_get_wtime() - time0
201
202 write(*,*)
203 write(*,*) 'End of Polygon example with 4 threads. Return code=',coi_error
204
205 If ( coi_error /= 0 ) then
206 call flog( "Quad threads: Errors encountered during solution", 1 )
207 elseif ( stacalls == 0 .or. solcalls == 0 ) then
208 call flog( "Quad threads: Status or Solution routine was not called", 1 )
209 elseif ( sstat /= 1 .or. mstat /= 2 ) then
210 call flog( "Quad threads: Solver and Model Status was not as expected (1,2)", 1 )
211 endif
212 iter4 = miter
213 xprim4 = xprim
214
215 same12 = .true.
216 same14 = .true.
217 same24 = .true.
218 Do i = 1, nv
219 If ( xprim2(i) /= xprim1(i) ) then
220 if ( same12 ) then
221 write(*,*) 'Solution 2 is different from Solution 1. First diff for index',i
222 endif
223 same12 = .false.
224 endif
225 If ( xprim4(i) /= xprim1(i) ) then
226 if ( same14 ) then
227 write(*,*) 'Solution 4 is different from Solution 1. First diff for index',i
228 endif
229 same14 = .false.
230 endif
231 If ( xprim4(i) /= xprim2(i) ) then
232 if ( same24 ) then
233 write(*,*) 'Solution 4 is different from Solution 2. First diff for index',i
234 endif
235 same24 = .false.
236 endif
237 Enddo
238
239 write(*,*)
240 write(*,"('Time for single thread',f10.3)") time1
241 write(*,"('Time for double thread',f10.3)") time2
242 write(*,"('Time for quad thread',f10.3)") time4
243 write(*,"('Speedup double ',f10.3)") time1/time2
244 write(*,"('Speedup quad ',f10.3)") time1/time4
245 if ( same12 .and. same14 .and. same24 ) then
246 write(*,*) 'All solutions are the same'
247 else
248 write(*,*) 'The solutions are NOT the same'
249! call flog( "Solutions are not the same", 1 )
250 Endif
251
252 if ( coi_free( cntvect ) /= 0 ) call flog( "Error while freeing control vector", 1 )
253
254 call flog( "Successful Solve", 0 )
255
256End Program polygon
257!
258! ============================================================================
259! Define information about the model:
260!
261
262!> Define information about the model
263!!
264!! @include{doc} readMatrix_params.dox
265Integer Function poly_readmatrix( lower, curr, upper, vsta, type, rhs, esta, &
266 colsta, rowno, value, nlflag, n, m, nz, &
267 usrmem )
268#if defined(itl)
269!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Poly_ReadMatrix
270#endif
271 Use pdata2
272 implicit none
273 integer, intent (in) :: n ! number of variables
274 integer, intent (in) :: m ! number of constraints
275 integer, intent (in) :: nz ! number of nonzeros
276 real*8, intent (in out), dimension(n) :: lower ! vector of lower bounds
277 real*8, intent (in out), dimension(n) :: curr ! vector of initial values
278 real*8, intent (in out), dimension(n) :: upper ! vector of upper bounds
279 integer, intent (in out), dimension(n) :: vsta ! vector of initial variable status
280 ! (not defined here)
281 integer, intent (out), dimension(m) :: type ! vector of equation types
282 integer, intent (in out), dimension(m) :: esta ! vector of initial equation status
283 ! (not defined here)
284 real*8, intent (in out), dimension(m) :: rhs ! vector of right hand sides
285 integer, intent (in out), dimension(n+1) :: colsta ! vector with start of column indices
286 integer, intent (out), dimension(nz) :: rowno ! vector of row numbers
287 integer, intent (in out), dimension(nz) :: nlflag ! vector of nonlinearity flags
288 real*8, intent (in out), dimension(nz) :: value ! vector of matrix values
289 real*8 usrmem(*) ! optional user memory
290
291 Integer :: i, j, k
292 real*8, parameter :: pi = 3.141592
293
294!
295! Information about Variables:
296! Default: Lower = -Inf, Curr = 0, and Upper = +inf.
297! Default: the status information in Vsta is not used.
298!
299! r.up(i) = 1;
300! theta.up(i) = pi;
301!
302! r.fx('i%nv%') = 0;
303! theta.fx('i%nv%') = pi;
304!
305! r.l(i) = 4*ord(i)*(card(i)+ 1- ord(i))/sqr(card(i)+1);
306! theta.l(i) = pi*ord(i)/card(i);
307 DO i = 1, np-1
308 lower(i) = 0.d0
309 lower(i+np) = 0.d0
310 upper(i) = 1.d0
311 upper(i+np) = pi
312 curr(i) = 4.d0*i*(np+1-i)/(1+np)**2
313 curr(i+np) = pi*i/np
314 enddo
315 lower(np) = 0.d0
316 upper(np) = 0.d0
317 curr(np) = 0.d0
318 lower(2*np) = pi
319 upper(2*np) = pi
320 curr(2*np) = pi
321!
322! Information about Constraints:
323! Default: Rhs = 0
324! Default: the status information in Esta and the function
325! value in FV are not used.
326! Default: Type: There is no default.
327! 0 = Equality,
328! 1 = Greater than or equal,
329! 2 = Less than or equal,
330! 3 = Non binding.
331!
332! ordered(i+1).. theta(i) =l= theta(i+1);
333!
334! distance(i,j)$(ord(j) > ord(i))..
335! sqr(r(i)) + sqr(r(j)) - 2*r(i)*r(j)*cos(theta(j)-theta(i)) =l= 1;
336!
337! obj.. polygon_area =e= 0.5*sum(j(i+1), r(i+1)*r(i)*sin(theta(i+1)-theta(i)));
338!
339 k = 0
340 DO i = 1, no
341 k = k + 1
342 Type(k) = 2 ! Ordered constraint Less than and RHS = 0
343 Enddo
344 rowmapk = 0 ! Initialize so we can recognize (i,j) combinations that are not equations
345 DO i = 1, np-1
346 DO j = i+1, np
347 k = k + 1
348 Type(k) = 2 ! Dist(i,j) Less than
349 rhs(k) = 1.d0
350 rowmapi(k-no) = i ! To be used for easier lookup
351 rowmapj(k-no) = j
352 rowmapk(i,j) = k
353 Enddo
354 Enddo
355 k = k + 1
356 type(k) = 3 ! Objective of type Non binding
357!
358! Information about the Jacobian. We have to define Rowno, Value,
359! Nlflag and Colsta.
360!
361! Colsta = Start of column indices (No Defaults):
362! Rowno = Row indices
363! Value = Value of derivative (by default only linear
364! derivatives are used)
365! Nlflag = 0 for linear and 1 for nonlinear derivative
366! (not needed for completely linear models)
367!
368! Indices
369!
370! r(1) r(2) r(3) r(4) .. Th(1) th(2) th(3) th(4)
371! 1: 1 -1
372! 2: 1 -1
373! 3: 1 -1
374! ..
375! NP: NL NL NL NL
376! NP+1: NL NL NL NL
377! NP+2 NL NL NL NL
378! ..
379! NL NL NL NL
380! NL NL NL NL
381! ..
382! Obj: NL NL NL NL NL NL
383!
384! Nonlinearity Structure: L = 0 are linear and NL = 1 are nonlinear
385! All nonzeros are nonlinear
386!
387 k = 0
388 DO i = 1, np ! r variables
389 colsta(i) = k+1
390 DO j = i+1, np
391 k = k + 1
392 rowno(k) = rowmapk(i,j) ! Distance - i index
393 nlflag(k) = 1
394 enddo
395 DO j = 1,np
396 IF ( rowmapk(j,i) > 0 ) Then
397 k = k + 1
398 rowno(k) = rowmapk(j,i) ! Distance - j index
399 nlflag(k) = 1
400 endif
401 enddo
402 k = k + 1
403 rowno(k) = no+nd+1 ! Objective constraint
404 nlflag(k) = 1
405 enddo
406 DO i = 1, np ! theta variables
407 colsta(np+i) = k+1
408 if ( i < np ) then
409 k = k + 1
410 rowno(k) = i ! ordered constraint index i
411 nlflag(k) = 0
412 value(k) = 1.d0
413 endif
414 if ( i > 1 ) then
415 k = k + 1
416 rowno(k) = i-1 ! ordered constraint index i+1
417 nlflag(k) = 0
418 value(k) = -1.d0
419 endif
420 DO j = i+1, np
421 k = k + 1
422 rowno(k) = rowmapk(i,j) ! Distance - i index
423 nlflag(k) = 1
424 enddo
425 DO j = 1,np
426 IF ( rowmapk(j,i) > 0 ) Then
427 k = k + 1
428 rowno(k) = rowmapk(j,i) ! Distance - j index
429 nlflag(k) = 1
430 endif
431 enddo
432 k = k + 1
433 rowno(k) = no+nd+1 ! Objective constraint
434 nlflag(k) = 1
435 enddo
436 colsta(nv+1) = k+1
437
438 poly_readmatrix = 0 ! Return value means OK
439
440end Function poly_readmatrix
441!
442!==========================================================================
443! Compute nonlinear terms and non-constant Jacobian elements
444!
445
446!> Compute nonlinear terms and non-constant Jacobian elements
447!!
448!! @include{doc} fdeval_params.dox
449Integer Function poly_fdeval( x, g, jac, rowno, jcnm, mode, ignerr, errcnt, &
450 n, nz, thread, usrmem )
451#if defined(itl)
452!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Poly_FDEval
453#endif
454 Use pdata2
455 implicit none
456 integer, intent (in) :: n ! number of variables
457 integer, intent (in) :: rowno ! number of the row to be evaluated
458 integer, intent (in) :: nz ! number of nonzeros in this row
459 real*8, intent (in), dimension(n) :: x ! vector of current solution values
460 real*8, intent (in out) :: g ! constraint value
461 real*8, intent (in out), dimension(n) :: jac ! vector of derivatives for current constraint
462 integer, intent (in), dimension(nz) :: jcnm ! list of variables that appear nonlinearly
463 ! in this row. Ffor information only.
464 integer, intent (in) :: mode ! evaluation mode: 1 = function value
465 ! 2 = derivatives, 3 = both
466 integer, intent (in) :: ignerr ! if 1 then errors can be ignored as long
467 ! as errcnt is incremented
468 integer, intent (in out) :: errcnt ! error counter to be incremented in case
469 ! of function evaluation errors.
470 integer, intent (in) :: thread
471 real*8 usrmem(*) ! optional user memory
472
473 Integer :: ne
474 Integer :: i, j
475 real*8 :: s, c
476
477 ne = n / 3
478!
479! Row NO+ND+1: the objective function is nonlinear
480!
481 if ( rowno == no+nd+1 ) then
482!
483! Mode = 1 or 3. Function value:
484! obj.. polygon_area =e= 0.5*sum(j(i+1), r(i+1)*r(i)*sin(theta(i+1)-theta(i)));
485!
486 if ( mode .eq. 1 .or. mode .eq. 3 ) then
487 g = 0.d0
488 do i = 1, np-1
489 g = g + x(i+1)*x(i)*sin(x(np+i+1)-x(np+i))
490 enddo
491 g = g * 0.5d0
492 endif
493!
494! Mode = 2 or 3: Derivative values: w.r.t. x: 2*(x(i)-x(j))*(-1/2)*dist(-3/2)
495!
496 if ( mode .eq. 2 .or. mode .eq. 3 ) then
497 do i = 1, nv
498 jac(i) = 0.d0
499 enddo
500 do i = 1, np-1
501! g = g + 0.5* x(i+1)*x(i)*Sin(x(NP+i+1)-x(NP+i))
502 s = 0.5d0*sin(x(np+i+1)-x(np+i))
503 jac(i+1) = jac(i+1) + x(i)*s
504 jac(i) = jac(i) + x(i+1)*s
505 c = 0.5d0*x(i+1)*x(i)*cos(x(np+i+1)-x(np+i))
506 jac(np+i+1) = jac(np+i+1) + c
507 jac(np+i) = jac(np+i) - c
508 enddo
509 endif
510!
511 Else if ( rowno > no ) then ! this is a distance constraint
512!
513! Mode = 1 or 3. Function value:
514! sqr(r(i)) + sqr(r(j)) - 2*r(i)*r(j)*cos(theta(j)-theta(i)) =l= 1;
515!
516 i = rowmapi(rowno-no)
517 j = rowmapj(rowno-no)
518 if ( mode .eq. 1 .or. mode .eq. 3 ) then
519 g = x(i)**2 + x(j)**2 -2.d0*x(i)*x(j)*cos(x(j+np)-x(i+np))
520 endif
521!
522! Mode = 2 or 3: Derivative values:
523!
524 if ( mode .eq. 2 .or. mode .eq. 3 ) then
525 c = cos(x(j+np)-x(i+np))
526 jac(i) = 2.d0*(x(i)-x(j)*c)
527 jac(j) = 2.d0*(x(j)-x(i)*c)
528 s = 2.d0*sin(x(j+np)-x(i+np))*x(i)*x(j)
529 jac(i+np) = -s
530 jac(j+np) = s
531 endif
532 Else
533 write(*,*) 'Fdeval called with rowno=',rowno
534 poly_fdeval = 1; return ! error in row number
535 endif
536 poly_fdeval = 0
537
538end Function poly_fdeval
539
540Integer Function poly_solution( XVAL, XMAR, XBAS, XSTA, YVAL, YMAR, YBAS, YSTA, N, M, USRMEM )
541#if defined(itl)
542!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Poly_Solution
543#endif
544!
545! Simple implementation in which we write the solution values to
546! the 'Documentation file' on unit 10.
547!
548 Use proginfop
549 Use pdata2
550 IMPLICIT NONE
551 INTEGER, Intent(IN) :: n, m
552 INTEGER, Intent(IN), Dimension(N) :: xbas, xsta
553 INTEGER, Intent(IN), Dimension(M) :: ybas, ysta
554 real*8, Intent(IN), Dimension(N) :: xval, xmar
555 real*8, Intent(IN), Dimension(M) :: yval, ymar
556 real*8, Intent(IN OUT) :: usrmem(*)
557
558 INTEGER i
559 CHARACTER*5, Parameter, Dimension(4) :: stat = (/ 'Lower','Upper','Basic','Super' /)
560
561 WRITE(10,"(/' Variable Solution value Reduced cost B-stat')")
562 DO i = 1, n
563 WRITE(10,"(1P,I7,E20.6,E16.6,4X,A5 )") i, xval(i), xmar(i), stat(1+xbas(i))
564 ENDDO
565
566 WRITE(10,"(/' Constrnt Activity level Marginal cost B-stat')")
567 DO i = 1, m
568 WRITE(10,"(1P,I7,E20.6,E16.6,4X,A5 )") i, yval(i), ymar(i), stat(1+ybas(i))
569 ENDDO
570 xprim(1:n) = xval(1:n)
571#if defined (wei)
572 Call flush(10)
573#endif
574
575 solcalls = solcalls + 1
576 poly_solution = 0
577
578END Function poly_solution
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_license(cntvect, licint1, licint2, licint3, licstring)
define the License Information.
Definition coistart.f90:680
integer function coidef_threadc(cntvect, threadc)
check for thread compatibility.
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
integer function poly_solution(xval, xmar, xbas, xsta, yval, ymar, ybas, ysta, n, m, usrmem)
program polygon
Main program. A simple setup and call of CONOPT.
integer, dimension(np, np) rowmapk
integer, parameter nd
integer, dimension(nd) rowmapi
real *8, dimension(nv) xprim
integer, dimension(nd) rowmapj
real *8, dimension(nv) xprim2
real *8, dimension(nv) xprim1
integer, parameter nv
integer, parameter no
real *8, dimension(nv) xprim4
integer, parameter np
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