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