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