CONOPT
Loading...
Searching...
No Matches
tutorial2r.f90
Go to the documentation of this file.
1!> @file tutorial2r.f90
2!! @ingroup FORT1THREAD_EXAMPLES
3!! @brief This is a revised version of Tutorial2 in which the production
4!! function is split into two equations using an intermediate variable.
5!!
6!! The purpose is to make the expressions simpler, in particular the
7!! second order terms.
8!! The equation
9!! Row 2: `Al*L**(-Rho) + Ak*K**(-Rho) + Ainp*Inp**(-Rho)) ** ( -1.d0/Rho ) - Out = 0`
10!! is replaced by
11!! Row 2 (new): `Al*L**(-Rho) + Ak*K**(-Rho) + Ainp*Inp**(-Rho)) = Int`
12!! Row 4 (new): `Int**( -1.d0/Rho ) = Out`
13!! where `Int` is a new variable.
14!!
15!! For more information about the individual callbacks, please have a look at the source code.
16
17#if defined(_WIN32) && !defined(_WIN64)
18#define dec_directives_win32
19#endif
20
21!> Main program. A simple setup and call of CONOPT
22!!
23Program tutorial2r
24
26 Use conopt
27 implicit None
28!
29! Declare the user callback routines as Integer, External:
30!
31 Integer, External :: tut_readmatrix ! Mandatory Matrix definition routine defined below
32 Integer, External :: tut_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 Integer, External :: tut_2dlagrstr ! Callback for second derivatives structure
39 Integer, External :: tut_2dlagrval ! Callback for second derivatives values
40#ifdef dec_directives_win32
41!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Tut_ReadMatrix
42!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Tut_FDEval
43!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Status
44!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Solution
45!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Message
46!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_ErrMsg
47!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Tut_2DLagrStr
48!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Tut_2DLagrVal
49#endif
50!
51! Control vector
52!
53 INTEGER, Dimension(:), Pointer :: cntvect
54 INTEGER :: coi_error
55
56 call startup
57!
58! Create and initialize a Control Vector
59!
60 coi_error = coi_create( cntvect )
61!
62! Tell CONOPT about the size of the model by populating the Control Vector:
63!
64 coi_error = max( coi_error, coidef_numvar( cntvect, 5 ) ) ! 5 variables
65 coi_error = max( coi_error, coidef_numcon( cntvect, 4 ) ) ! 4 constraints
66 coi_error = max( coi_error, coidef_numnz( cntvect, 11 ) ) ! 11 nonzeros in the Jacobian
67 coi_error = max( coi_error, coidef_numnlnz( cntvect, 5 ) ) ! 5 of which are nonlinear
68 coi_error = max( coi_error, coidef_numhess( cntvect, 4 ) ) ! 4 elements in the Hessian of the Lagrangian
69 coi_error = max( coi_error, coidef_optdir( cntvect, 1 ) ) ! Maximize
70 coi_error = max( coi_error, coidef_objcon( cntvect, 1 ) ) ! Objective is constraint 1
71 coi_error = max( coi_error, coidef_debugfv( cntvect, 0 ) ) ! Debug Fdeval on or off
72 coi_error = max( coi_error, coidef_debug2d( cntvect, 0 ) ) ! Debug 2nd derivatives on or off
73 coi_error = max( coi_error, coidef_optfile( cntvect, 'tutorial.opt' ) )
74!
75! Tell CONOPT about the callback routines:
76!
77 coi_error = max( coi_error, coidef_readmatrix( cntvect, tut_readmatrix ) )
78 coi_error = max( coi_error, coidef_fdeval( cntvect, tut_fdeval ) )
79 coi_error = max( coi_error, coidef_status( cntvect, std_status ) )
80 coi_error = max( coi_error, coidef_solution( cntvect, std_solution ) )
81 coi_error = max( coi_error, coidef_message( cntvect, std_message ) )
82 coi_error = max( coi_error, coidef_errmsg( cntvect, std_errmsg ) )
83 coi_error = max( coi_error, coidef_2dlagrstr( cntvect, tut_2dlagrstr ) )
84 coi_error = max( coi_error, coidef_2dlagrval( cntvect, tut_2dlagrval ) )
85
86#if defined(CONOPT_LICENSE_INT_1) && defined(CONOPT_LICENSE_INT_2) && defined(CONOPT_LICENSE_INT_3) && defined(CONOPT_LICENSE_TEXT)
87 coi_error = max( coi_error, coidef_license( cntvect, conopt_license_int_1, conopt_license_int_2, conopt_license_int_3, conopt_license_text) )
88#endif
89
90 If ( coi_error .ne. 0 ) THEN
91 write(*,*)
92 write(*,*) '**** Fatal Error while loading CONOPT Callback routines.'
93 write(*,*)
94 call flog( "Skipping Solve due to setup errors", 1 )
95 ENDIF
96!
97! Start CONOPT:
98!
99 coi_error = coi_solve( cntvect )
100
101 write(*,*)
102 write(*,*) 'End of Tutorial2r example. Return code=',coi_error
103
104 If ( coi_error /= 0 ) then
105 call flog( "Errors encountered during solution", 1 )
106 elseif ( stacalls == 0 .or. solcalls == 0 ) then
107 call flog( "Status or Solution routine was not called", 1 )
108 elseif ( sstat /= 1 .or. mstat /= 2 ) then
109 call flog( "Solver and Model Status was not as expected (1,2)", 1 )
110 elseif ( abs( obj-0.572943d0 ) > 0.000001d0 ) then
111 call flog( "Incorrect objective returned", 1 )
112 endif
113
114 if ( coi_free(cntvect) /= 0 ) call flog( "Error while freeing control vector",1)
115
116 call flog( "Successful Solve", 0 )
117
118End Program tutorial2r
119
120!> Define information about the model
121!!
122!! @include{doc} readMatrix_params.dox
123Integer Function tut_readmatrix( lower, curr, upper, vsta, type, rhs, esta, &
124 colsta, rowno, value, nlflag, n, m, nz, &
125 usrmem )
126#ifdef dec_directives_win32
127!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Tut_ReadMatrix
128#endif
129 implicit none
130 integer, intent (in) :: n ! number of variables
131 integer, intent (in) :: m ! number of constraints
132 integer, intent (in) :: nz ! number of nonzeros
133 real*8, intent (in out), dimension(n) :: lower ! vector of lower bounds
134 real*8, intent (in out), dimension(n) :: curr ! vector of initial values
135 real*8, intent (in out), dimension(n) :: upper ! vector of upper bounds
136 integer, intent (in out), dimension(n) :: vsta ! vector of initial variable status
137 ! (not defined here)
138 integer, intent (out), dimension(m) :: type ! vector of equation types
139 integer, intent (in out), dimension(m) :: esta ! vector of initial equation status
140 ! (not defined here)
141 real*8, intent (in out), dimension(m) :: rhs ! vector of right hand sides
142 integer, intent (in out), dimension(n+1) :: colsta ! vector with start of column indices
143 integer, intent (out), dimension(nz) :: rowno ! vector of row numbers
144 integer, intent (in out), dimension(nz) :: nlflag ! vector of nonlinearity flags
145 real*8, intent (in out), dimension(nz) :: value ! vector of matrix values
146 real*8 usrmem(*) ! optional user memory
147!
148! Parameters needed for new initial value -- could be moved to a
149! Module
150!
151 real*8, parameter :: al = 0.16d0
152 real*8, parameter :: ak = 2.0d0
153 real*8, parameter :: ainp = 0.16d0
154 real*8, parameter :: rho = 1.0d0
155 real*8, parameter :: k = 4.0d0
156!
157! Information about Variables:
158! Default: Lower = -Inf, Curr = 0, and Upper = +inf.
159! Default: the status information in Vsta is not used.
160!
161! Lower bound on L = X(1) = 0.1 and initial value = 0.5:
162!
163 lower(1) = 0.1d0
164 curr(1) = 0.5d0
165!
166! Lower bound on INP = X(2) = 0.1 and initial value = 0.5:
167!
168 lower(2) = 0.1d0
169 curr(2) = 0.5d0
170!
171! Lower bound on OUT = X(3) and P = X(4) are both 0 and the
172! default initial value of 0 is used:
173!
174 lower(3) = 0.d0
175 lower(4) = 0.d0
176!
177! The added variable, in the following called Int, is given a lower
178! bound of 0.01 and an initial value consistent with the defining
179! equation, i.e. int = (Al*L**(-Rho) + Ak*K**(-Rho) + Ainp*Inp**(-Rho))
180!
181 lower(5) = 0.01d0
182 curr(5) = (al*curr(1)**(-rho) + ak*k**(-rho) + ainp*curr(2)**(-rho))
183!
184! Information about Constraints:
185! Default: Rhs = 0
186! Default: the status information in Esta and the function
187! value in FV are not used.
188! Default: Type: There is no default.
189! 0 = Equality,
190! 1 = Greater than or equal,
191! 2 = Less than or equal,
192! 3 = Non binding.
193!
194! Constraint 1 (Objective)
195! Rhs = -0.1 and type Non binding
196!
197 rhs(1) = -0.1d0
198 type(1) = 3
199!
200! Constraint 2 (Production Function)
201! Rhs = 0 and type Equality
202!
203 type(2) = 0
204!
205! Constraint 3 (Price equation)
206! Rhs = 4.0 and type Equality
207!
208 rhs(3) = 4.d0
209 type(3) = 0
210!
211! Constraint 4 (additional equation)
212! Rhs = 0.0 and type Equality
213!
214 type(4) = 0
215!
216! Information about the Jacobian. CONOPT expects a columnwise
217! representation in Rowno, Value, Nlflag and Colsta.
218!
219! Colsta = Start of column indices (No Defaults):
220! Rowno = Row indices
221! Value = Value of derivative (by default only linear
222! derivatives are used)
223! Nlflag = 0 for linear and 1 for nonlinear derivative
224! (not needed for completely linear models)
225!
226! Indices
227! L INP OUT P Int
228! x(1) x(2) x(3) x(4) x(5)
229! 1: 1 3 5 8
230! 2: 2 4 10
231! 3: 6 9
232! 4 7 11
233!
234 colsta(1) = 1
235 colsta(2) = 3
236 colsta(3) = 5
237 colsta(4) = 8
238 colsta(5) = 10
239 colsta(6) = 12
240 rowno(1) = 1
241 rowno(2) = 2
242 rowno(3) = 1
243 rowno(4) = 2
244 rowno(5) = 1
245 rowno(6) = 3
246 rowno(7) = 4
247 rowno(8) = 1
248 rowno(9) = 3
249 rowno(10) = 2
250 rowno(11) = 4
251!
252! Nonlinearity Structure: L = 0 are linear and NL = 1 are nonlinear
253! x(1) x(2) x(3) x(4) x(5)
254! 1: L L NL NL
255! 2: NL NL L
256! 3: L L
257! 4: L NL
258!
259 nlflag(1) = 0
260 nlflag(2) = 1
261 nlflag(3) = 0
262 nlflag(4) = 1
263 nlflag(5) = 1
264 nlflag(6) = 0
265 nlflag(7) = 0
266 nlflag(8) = 1
267 nlflag(9) = 0
268 nlflag(10)= 0
269 nlflag(11)= 1
270!
271! Value (Linear only)
272! x(1) x(2) x(3) x(4) x(5)
273! 1: -1 -1 NL NL
274! 2: NL NL -1
275! 3: 1 2
276! 4: -1 NL
277!
278 value(1) = -1.d0
279 value(3) = -1.d0
280 value(6) = 1.d0
281 value(7) = -1.d0
282 value(9) = 2.d0
283 value(10) = -1.d0
284
285 tut_readmatrix = 0 ! Return value means OK
286
287end Function tut_readmatrix
288
289!> Compute nonlinear terms and non-constant Jacobian elements
290!!
291!! @include{doc} fdeval_params.dox
292Integer Function tut_fdeval( x, g, jac, rowno, jcnm, mode, ignerr, errcnt, &
293 n, nz, thread, usrmem )
294#ifdef dec_directives_win32
295!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Tut_FDEval
296#endif
297 implicit none
298 integer, intent (in) :: n ! number of variables
299 integer, intent (in) :: rowno ! number of the row to be evaluated
300 integer, intent (in) :: nz ! number of nonzeros in this row
301 real*8, intent (in), dimension(n) :: x ! vector of current solution values
302 real*8, intent (in out) :: g ! constraint value
303 real*8, intent (in out), dimension(n) :: jac ! vector of derivatives for current constraint
304 integer, intent (in), dimension(nz) :: jcnm ! list of variables that appear nonlinearly
305 ! in this row. Ffor information only.
306 integer, intent (in) :: mode ! evaluation mode: 1 = function value
307 ! 2 = derivatives, 3 = both
308 integer, intent (in) :: ignerr ! if 1 then errors can be ignored as long
309 ! as errcnt is incremented
310 integer, intent (in out) :: errcnt ! error counter to be incremented in case
311 ! of function evaluation errors.
312 integer, intent (in) :: thread
313 real*8 usrmem(*) ! optional user memory
314!
315! Declare local copies of the optimization variables. This is
316! just for convenience to make the expressions easier to read.
317!
318 real*8 :: l, inp, out, p, int
319!
320! Declare parameters and their data values.
321!
322 real*8, parameter :: w = 1.0d0
323 real*8, parameter :: l0 = 0.1d0
324 real*8, parameter :: pinp = 1.0d0
325 real*8, parameter :: al = 0.16d0
326 real*8, parameter :: ak = 2.0d0
327 real*8, parameter :: ainp = 0.16d0
328 real*8, parameter :: rho = 1.0d0
329 real*8, parameter :: k = 4.0d0
330!
331! Move the optimization variables from the X vector to a set
332! of local variables with the same names as the variables in
333! the model description. This is not necessary, but it should make
334! the equations easier to recognize.
335!
336 l = x(1)
337 inp = x(2)
338 out = x(3)
339 p = x(4)
340 int = x(5)
341!
342! Row 1: the objective function is nonlinear
343!
344 if ( rowno .eq. 1 ) then
345!
346! Mode = 1 or 3. Function value: G = P * Out
347!
348 if ( mode .eq. 1 .or. mode .eq. 3 ) then
349 g = p * out
350 endif
351!
352! Mode = 2 or 3: Derivative values:
353!
354 if ( mode .eq. 2 .or. mode .eq. 3 ) then
355 jac(3) = p ! derivative w.r.t. Out = X(3)
356 jac(4) = out ! derivative w.r.t. P = X(4)
357 endif
358!
359! Row 2: The production function is nonlinear:
360! Al*L**(-Rho) + Ak*K**(-Rho) + Ainp*Inp**(-Rho))
361!
362 elseif ( rowno .eq. 2 ) then
363!
364! Mode = 1 or 3: Function value
365!
366 if ( mode .eq. 1 .or. mode .eq. 3 ) then
367 g = al*l**(-rho) + ak*k**(-rho) + ainp*inp**(-rho)
368 endif
369!
370! Mode = 2 or 3: Derivatives
371!
372 if ( mode .eq. 2 .or. mode .eq. 3 ) then
373 jac(1) = -rho * al * l ** (-rho-1.d0) ! derivative w.r.t. L = X(1)
374 jac(2) = -rho * ainp * inp ** (-rho-1.d0) ! derivative w.r.t. Inp = X(2)
375 endif
376!
377! Row = 4: Out = Int**(-1/Rho)
379 elseif ( rowno .eq. 4 ) then
380 if ( mode .eq. 1 .or. mode .eq. 3 ) then
381 g = int**(-1.d0/rho)
382 endif
383 if ( mode .eq. 2 .or. mode .eq. 3 ) then
384 jac(5) = (-1.d0/rho)*int**(-1.d0/rho-1.d0) ! derivative w.r.t. Int = X(5)
385 endif
386
387 endif
388 tut_fdeval = 0
389
390end Function tut_fdeval
391
392!> Specify the structure of the Lagrangian of the Hessian
393!!
394!! @include{doc} 2DLagrStr_params.dox
395Integer Function tut_2dlagrstr( HSRW, HSCL, NODRV, N, M, NHESS, UsrMem )
396#ifdef dec_directives_win32
397!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Tut_2DLagrStr
398#endif
399 Implicit None
400
401 Integer, Intent (IN) :: n, m, nhess
402 Integer, Intent (IN OUT) :: nodrv
403 Integer, Dimension(Nhess), Intent (Out) :: hsrw, hscl
404 real*8, Intent(IN OUT) :: usrmem(*)
405!
406! There are two nonlinear constraints, P * Out in row 1 and Sdef and
407! Al*L**(-Rho) + Ak*K**(-Rho) + Ainp*Inp**(-Rho)) ** ( -1.d0/Rho ) in
408! row 2 and they have the following Hessians:
409!
410! Row1: P * Out
411! Hessian (diagonal and lower triangle): A total of 1 nonzero element
412! Out P
413! 3 4
414! 3
415! 4 1
416!
417! Row 2: Al*L**(-Rho) + Ak*K**(-Rho) + Ainp*Inp**(-Rho))
418! with L and Inp being the variables.
419! Hessian (diagonal and lower triangle): A total of 3 nonzero element
420! L Inp
421! 1 2
422! 1 d1
423! 2 d2
424!
425! where d1 = (-Rho)*(-Rho-1)*Al*L**(-Rho-2)
426! and d2 = (-Rho)*(-Rho-1)*Ainp*Inp**(-Rho-2)
427!
428! Row 4: Int**(-1/Rho)
429! Hessian (diagonal and lower triangle): A total of 1 nonzero element
430! Int
431! 5
432! 5 Int d3
433!
434! where d3 = (-1/Rho)*(-1/Rho-1)*Int**(-1/Rho-2)
435!
436! The three Hessian matrices do not overlap so the total number of
437! nonzeros is 4. The structure using numerical indices for rows and
438! columns and running indices for the Hessian element is:
439!
440! 1 2 3 4 5
441! 1 1
442! 2 2
443! 3
444! 4 3
445! 5 4
446!
447! and the same with names of variables and numbers of contributing
448! equations:
449!
450! L Inp Out P Int
451! L 2
452! Inp 2
453! Out
454! P 1
455! Int 4
456!
457!
458! Define structure of Hessian
459!
460 hsrw(1) = 1; hscl(1) = 1
461 hsrw(2) = 2; hscl(2) = 2
462 hsrw(3) = 4; hscl(3) = 3
463 hsrw(4) = 5; hscl(4) = 5
464
465 tut_2dlagrstr = 0
466
467End Function tut_2dlagrstr
468
469!> Compute the Lagrangian of the Hessian
470!!
471!! @include{doc} 2DLagrVal_params.dox
472Integer Function tut_2dlagrval( X, U, HSRW, HSCL, HSVL, NODRV, N, M, NHESS, UsrMem )
473#ifdef dec_directives_win32
474!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Tut_2DLagrVal
475#endif
476 Implicit None
477
478 Integer, Intent (IN) :: n, m, nhess
479 Integer, Intent (IN OUT) :: nodrv
480 real*8, Dimension(N), Intent (IN) :: x
481 real*8, Dimension(M), Intent (IN) :: u
482 Integer, Dimension(NHess), Intent (IN) :: hsrw, hscl
483 real*8, Dimension(NHess), Intent (Out) :: hsvl
484 real*8, Intent(IN OUT) :: usrmem(*)
485!
486! Declare local copies of the optimization variables. This is
487! just for convenience to make the expressions easier to read.
488!
489 real*8 :: l, inp
490!
491! Declare parameters and their data values.
492!
493 real*8, parameter :: w = 1.0d0
494 real*8, parameter :: l0 = 0.1d0
495 real*8, parameter :: pinp = 1.0d0
496 real*8, parameter :: al = 0.16d0
497 real*8, parameter :: ak = 2.0d0
498 real*8, parameter :: ainp = 0.16d0
499 real*8, parameter :: rho = 1.0d0
500 real*8, parameter :: k = 4.0d0
501!
502! There are two nonlinear constraints, P * Out in row 1 and Sdef and
503! Al*L**(-Rho) + Ak*K**(-Rho) + Ainp*Inp**(-Rho)) ** ( -1.d0/Rho ) in
504! row 2 and they have the following Hessians:
505!
506! Row1: P * Out
507! Hessian (diagonal and lower triangle): A total of 1 nonzero element
508! Out P
509! 3 4
510! 3
511! 4 1
512!
513! Row 2: Al*L**(-Rho) + Ak*K**(-Rho) + Ainp*Inp**(-Rho))
514! with L and Inp being the variables.
515! Hessian (diagonal and lower triangle): A total of 3 nonzero element
516! L Inp
517! 1 2
518! 1 d1
519! 2 d2
520!
521! where d1 = (-Rho)*(-Rho-1)*Al*L**(-Rho-2)
522! and d2 = (-Rho)*(-Rho-1)*Ainp*Inp**(-Rho-2)
523!
524! Row 4: Int**(-1/Rho)
525! Hessian (diagonal and lower triangle): A total of 1 nonzero element
526! Int
527! 5
528! 5 Int d3
529!
530! where d3 = (-1/Rho)*(-1/Rho-1)*Int**(-1/Rho-2)
531!
532! The three Hessian matrices do not overlap so the total number of
533! nonzeros is 4. The structure using numerical indices for rows and
534! columns and running indices for the Hessian element is:
535!
536! 1 2 3 4 5
537! 1 1
538! 2 2
539! 3
540! 4 3
541! 5 4
542!
543! and the same with names of variables and numbers of contributing
544! equations:
545!
546! L Inp Out P Int
547! L 2
548! Inp 2
549! Out
550! P 1
551! Int 4
552!
553!
554! Normal Evaluation mode
555!
556 hsvl(3) = u(1) ! Second derivative of constraint 1 is 1
557 if ( u(2) .ne. 0.d0 ) Then
558 l = x(1)
559 inp = x(2)
560 hsvl(1) = (-rho) * (-rho-1.d0) * al * l**(-rho-2.d0) * u(2)
561 hsvl(2) = (-rho) * (-rho-1.d0) * ainp * inp**(-rho-2.d0) * u(2)
562 endif
563 if ( u(4) .ne. 0.d0 ) then
564 hsvl(4) = (-1.d0/rho)*(-1.d0/rho-1.d0)*x(5)**(-1.d0/rho-2.d0)
565 endif
566
567 tut_2dlagrval = 0
568
569End Function tut_2dlagrval
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_2dlagrstr(cntvect, coi_2dlagrstr)
define callback routine for providing the structure of the second derivatives of the Lagrangian.
Definition conopt.f90:1547
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_2dlagrval(cntvect, coi_2dlagrval)
define callback routine for computing the values of the second derivatives of the Lagrangian.
Definition conopt.f90:1573
integer(c_int) function coidef_license(cntvect, licint1, licint2, licint3, licstring)
define the License Information.
Definition conopt.f90:293
integer(c_int) function coidef_debugfv(cntvect, debugfv)
turn Debugging of FDEval on and off.
Definition conopt.f90:387
integer(c_int) function coidef_debug2d(cntvect, debug2d)
turn debugging of 2nd derivatives on and off.
Definition conopt.f90:534
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_numhess(cntvect, numhess)
defines the Number of Hessian Nonzeros.
Definition conopt.f90:188
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
real *8 obj
Definition comdecl.f90:16
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 tut_2dlagrval(x, u, hsrw, hscl, hsvl, nodrv, n, m, nhess, usrmem)
Compute the Lagrangian of the Hessian.
integer function tut_2dlagrstr(hsrw, hscl, nodrv, n, m, nhess, usrmem)
Specify the structure of the Lagrangian of the Hessian.
program tutorial2r
Main program. A simple setup and call of CONOPT.
integer function tut_readmatrix(lower, curr, upper, vsta, type, rhs, esta, colsta, rowno, value, nlflag, n, m, nz, usrmem)
Define information about the model.
Definition tutorial.f90:104
integer function tut_fdeval(x, g, jac, rowno, jcnm, mode, ignerr, errcnt, n, nz, thread, usrmem)
Compute nonlinear terms and non-constant Jacobian elements.
Definition tutorial.f90:237