CONOPT
Loading...
Searching...
No Matches
leastsq6.f90
Go to the documentation of this file.
1!> @file leastsq6.f90
2!! @ingroup FORT1THREAD_EXAMPLES
3!!
4!! Currently, the model crashes in the inversion -- have it fixed later.
5!!
6!! The model is similar to leastsq, but we test that the duals have
7!! the right value, run again without a post triangle and check the
8!! duals and we also compare duals with the change in objective from a
9!! perturbation.
10!!
11!! We solve the following nonlinear least squares model:
12!!
13!! \f[
14!! \min \sum_i res_{i}^2 \\
15!!
16!! \sum_j ( a_{ij}x_j + b_{ij}x_j^2 ) + res_i = obs_i
17!! \f]
18!!
19!! where \f$a\f$, \f$b\f$, and \f$obs\f$ are known data, and \f$res\f$ and \f$x\f$ are the
20!! variables of the model.
21!!
22!!
23!! Casenum = 1 standard,
24!! Casenum = 2 with wide bounds on Res
25!! Casenum = 3 standard model with obs(1) changed by +1.e-3
26!!
27!! For more information about the individual callbacks, please have a look at the source code.
28
29#if defined(_WIN32) && !defined(_WIN64)
30#define dec_directives_win32
31#endif
32
33REAL FUNCTION rndx( )
34!
35! Defines a pseudo random number between 0 and 1
36!
37 IMPLICIT NONE
38
39 Integer, save :: seed = 12359
40
41 seed = mod(seed*1027+25,1048576)
42 rndx = float(seed)/float(1048576)
43
44END FUNCTION rndx
45
46subroutine defdata
47!
48! Define values for A, B, and Obs
49!
50 Use lsq_70
51 IMPLICIT NONE
52
53 Integer :: i, j
54 real*8, Parameter :: xtarg = -1.0
55 real*8, Parameter :: noise = 1.0
56 Real, External :: Rndx
57 real*8 :: o
58
59 do i = 1, nobs
60 o = 0.d0
61 do j = 1, dimx
62 a(i,j) = rndx()
63 b(i,j) = rndx()
64 o = o + a(i,j) * xtarg + b(i,j) * xtarg**2
65 enddo
66 obs(i) = o + noise * rndx()
67 enddo
68
69end subroutine defdata
70!
71! Main program.
72!
73!> Main program. A simple setup and call of CONOPT
74!!
75Program leastsquare
76
78 Use conopt
79 Use lsq_70
80 Use casedata_num
81 implicit None
82!
83! Declare the user callback routines as Integer, External:
84!
85 Integer, External :: lsq_readmatrix ! Mandatory Matrix definition routine defined below
86 Integer, External :: lsq_fdeval ! Function and Derivative evaluation routine
87 ! needed a nonlinear model.
88 Integer, External :: std_status ! Standard callback for displaying solution status
89 Integer, External :: std_solution ! Standard callback for displaying solution values
90 Integer, External :: std_message ! Standard callback for managing messages
91 Integer, External :: std_errmsg ! Standard callback for managing error messages
92 Integer, External :: lsq_2dlagrstr ! Callback for second derivatives structure
93 Integer, External :: lsq_2dlagrval ! Callback for second derivatives values
94 Integer, External :: lsq_option ! Option defining callback routine
95#ifdef dec_directives_win32
96!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_ReadMatrix
97!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_FDEval
98!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Status
99!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Solution
100!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Message
101!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_ErrMsg
102!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_2DLagrStr
103!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_2DLagrVal
104!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_Option
105#endif
106!
107! Control vector
108!
109 INTEGER, Dimension(:), Pointer :: cntvect
110 INTEGER :: coi_error
111!
112 Integer :: i
113 Character(Len=132) :: text
114 real*8 :: obj_run1, obj_run3
115
116 call startup
117!
118! Define data
119!
120 Call defdata
121!
122! Create and initialize a Control Vector
123!
124 coi_error = coi_create( cntvect )
125!
126! Tell CONOPT about the size of the model by populating the Control Vector:
127!
128 coi_error = max( coi_error, coidef_numvar( cntvect, nobs + dimx ) )
129 coi_error = max( coi_error, coidef_numcon( cntvect, nobs + 1 ) )
130 coi_error = max( coi_error, coidef_numnz( cntvect, nobs * (dimx+2) ) )
131 coi_error = max( coi_error, coidef_numnlnz( cntvect, nobs * (dimx+1) ) )
132 coi_error = max( coi_error, coidef_numhess( cntvect, nobs + dimx ) )
133 coi_error = max( coi_error, coidef_optdir( cntvect, -1 ) ) ! Minimize
134 coi_error = max( coi_error, coidef_objcon( cntvect, nobs + 1 ) ) ! Objective is last constraint
135 coi_error = max( coi_error, coidef_optfile( cntvect, 'leastsq6.opt' ) )
136!
137! Tell CONOPT about the callback routines:
138!
139 coi_error = max( coi_error, coidef_readmatrix( cntvect, lsq_readmatrix ) )
140 coi_error = max( coi_error, coidef_fdeval( cntvect, lsq_fdeval ) )
141 coi_error = max( coi_error, coidef_status( cntvect, std_status ) )
142 coi_error = max( coi_error, coidef_solution( cntvect, std_solution ) )
143 coi_error = max( coi_error, coidef_message( cntvect, std_message ) )
144 coi_error = max( coi_error, coidef_errmsg( cntvect, std_errmsg ) )
145 coi_error = max( coi_error, coidef_2dlagrstr( cntvect, lsq_2dlagrstr ) )
146 coi_error = max( coi_error, coidef_2dlagrval( cntvect, lsq_2dlagrval ) )
147 coi_error = max( coi_error, coidef_option( cntvect, lsq_option ) )
148
149#if defined(CONOPT_LICENSE_INT_1) && defined(CONOPT_LICENSE_INT_2) && defined(CONOPT_LICENSE_INT_3) && defined(CONOPT_LICENSE_TEXT)
150 coi_error = max( coi_error, coidef_license( cntvect, conopt_license_int_1, conopt_license_int_2, conopt_license_int_3, conopt_license_text) )
151#endif
152
153 If ( coi_error .ne. 0 ) THEN
154 write(*,*)
155 write(*,*) '**** Fatal Error while loading CONOPT Callback routines.'
156 write(*,*)
157 call flog( "Skipping Solve due to setup errors", 1 )
158 ENDIF
159!
160! Save the solution so we can check the duals:
161!
162 do_allocate = .true.
163!
164! Start CONOPT:
165!
166 casenum = 1
167 coi_error = coi_solve( cntvect )
168
169 If ( coi_error /= 0 ) then
170 call flog( "Run 1: Errors encountered during solution", 1 )
171 elseif ( stacalls == 0 .or. solcalls == 0 ) then
172 call flog( "Run 1: Status or Solution routine was not called", 1 )
173 elseif ( .not. ( sstat == 1 .and. mstat == 2 ) ) then
174 call flog( "Run 1: Solver or Model status was not as expected (1,2)", 1 )
175 Else ! Objective is not know and therefore not tested
176 do i = 1, nobs
177 if ( abs( udual(i) - 2*xprim(dimx+i) ) > 1.d-6 ) then
178 write(10,*) 'Run 1: Dual and Primal did not match in constraint',i,' P=',xprim(dimx+i),' D=',udual(i),' Error=', &
179 abs( udual(i) - 2*xprim(dimx+i) )
180 write(text,*) 'Run 1: Dual and Primal did not match in constraint',i,' P=',xprim(dimx+i),' D=',udual(i)
181 call flog( text, 1 )
182 endif
183 enddo
184 Call checkdual( 'Run 1: Leastsq6', minimize )
185 endif
186 obj_run1 = obj
187 obj_run3 = obj + 1.0d-3*udual(1) ! expected new objective in Run 3
188!
189! Repeat but with bounds on the residuals. They cannot be moved into
190! the post-triangle any more.
191!
192 casenum = 2
193 coi_error = coi_solve( cntvect )
194
195 If ( coi_error /= 0 ) then
196 call flog( "Run 2: Errors encountered during solution", 1 )
197 elseif ( stacalls == 0 .or. solcalls == 0 ) then
198 call flog( "Run 2: Status or Solution routine was not called", 1 )
199 elseif ( .not. ( sstat == 1 .and. mstat == 2 ) ) then
200 call flog( "Run 2: Solver or Model status was not as expected (1,2)", 1 )
201 elseif ( abs( obj - obj_run1 ) > 1.d-7 ) then
202 call flog( "Run 2: Incorrect objective returned", 1 )
203 Else
204 do i = 1, nobs
205 if ( abs( udual(i) - 2*xprim(dimx+i) ) > 1.d-6 ) then
206 write(10,*) 'Run 2: Dual and Primal did not match in constraint',i,' P=',xprim(dimx+i),' D=',udual(i),' Error=', &
207 abs( udual(i) - 2*xprim(dimx+i) )
208 write(text,*) 'Run 2: Dual and Primal did not match in constraint',i,' P=',xprim(dimx+i),' D=',udual(i)
209 call flog( text, 1 )
210 endif
211 enddo
212 Call checkdual( 'Run 2: Leastsq6', minimize )
213 endif
214!
215! Repeat but with Obs(1) changed by 1.d-3. The objective should change
216! by 1.d-3*Udual(1).
217!
218 casenum = 3
219 obs(1) = obs(1) + 1.0d-3
220 write(10,*) 'Actual Objective for run 2=',obj
221 write(10,*) 'Expected Objective for run 3=',obj_run3
222 coi_error = coi_solve( cntvect )
223 write(10,*) 'Actual Objective for run 3=',obj
224 write(10,*) 'Expected Objective for run 3=',obj_run3
225 write(10,*) 'Diffrens Objective for run 3=', abs(obj_run3-obj)
226
227 If ( coi_error /= 0 ) then
228 call flog( "Run 3: Errors encountered during solution", 1 )
229 elseif ( stacalls == 0 .or. solcalls == 0 ) then
230 call flog( "Run 3: Status or Solution routine was not called", 1 )
231 elseif ( .not. ( sstat == 1 .and. mstat == 2 ) ) then
232 call flog( "Run 3: Solver or Model status was not as expected (1,2)", 1 )
233 elseif ( abs( obj - obj_run3 ) > 1.d-6 ) then
234 call flog( "Run 3: Incorrect objective returned in run3 ", 1 )
235 else
236 Call checkdual( 'Run 3: Leastsq6', minimize )
237 endif
238
239 write(*,*)
240 write(*,*) 'End of Least Square example 6. Return code=',coi_error
241
242 if ( coi_free(cntvect) /= 0 ) call flog( "Error while freeing control vector",1)
244 call flog( "Successful Solve", 0 )
245
246End Program leastsquare
247!
248! ============================================================================
249! Define information about the model:
250!
251
252!> Define information about the model
253!!
254!! @include{doc} readMatrix_params.dox
255Integer Function lsq_readmatrix( lower, curr, upper, vsta, type, rhs, esta, &
256 colsta, rowno, value, nlflag, n, m, nz, &
257 usrmem )
258#ifdef dec_directives_win32
259!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_ReadMatrix
260#endif
261 Use lsq_70
262 Use casedata_num
263 implicit none
264 integer, intent (in) :: n ! number of variables
265 integer, intent (in) :: m ! number of constraints
266 integer, intent (in) :: nz ! number of nonzeros
267 real*8, intent (in out), dimension(n) :: lower ! vector of lower bounds
268 real*8, intent (in out), dimension(n) :: curr ! vector of initial values
269 real*8, intent (in out), dimension(n) :: upper ! vector of upper bounds
270 integer, intent (in out), dimension(n) :: vsta ! vector of initial variable status
271 ! (not defined here)
272 integer, intent (out), dimension(m) :: type ! vector of equation types
273 integer, intent (in out), dimension(m) :: esta ! vector of initial equation status
274 ! (not defined here)
275 real*8, intent (in out), dimension(m) :: rhs ! vector of right hand sides
276 integer, intent (in out), dimension(n+1) :: colsta ! vector with start of column indices
277 integer, intent (out), dimension(nz) :: rowno ! vector of row numbers
278 integer, intent (in out), dimension(nz) :: nlflag ! vector of nonlinearity flags
279 real*8, intent (in out), dimension(nz) :: value ! vector of matrix values
280 real*8 usrmem(*) ! optional user memory
281
282 Integer :: i, j, k
283!
284! Information about Variables:
285! Default: Lower = -Inf, Curr = 0, and Upper = +inf.
286! Default: the status information in Vsta is not used.
287!
288 do i = 1, dimx
289 curr(i) = -0.8d0
290 enddo
291 if ( casenum == 2 ) Then ! Add large bound on Res that will prevent the
292 do i = 1, nobs ! post triangle from being identified.
293 lower(dimx+i) = -1000.0d0
294 upper(dimx+i) = +1000.0d0
295 enddo
296 endif
297!
298! Information about Constraints:
299! Default: Rhs = 0
300! Default: the status information in Esta and the function
301! value in FV are not used.
302! Default: Type: There is no default.
303! 0 = Equality,
304! 1 = Greater than or equal,
305! 2 = Less than or equal,
306! 3 = Non binding.
307!
308! Constraints 1 to Nobs:
309! Rhs = Obs(i) and type Equality
310!
311 do i = 1, nobs
312 rhs(i) = obs(i)
313 type(i) = 0
314 enddo
315!
316! Constraint Nobs + 1 (Objective)
317! Rhs = 0 and type Non binding
318!
319 type(nobs+1) = 3
320!
321! Information about the Jacobian. CONOPT expects a columnwise
322! representation in Rowno, Value, Nlflag and Colsta.
323!
324! Colsta = Start of column indices (No Defaults):
325! Rowno = Row indices
326! Value = Value of derivative (by default only linear
327! derivatives are used)
328! Nlflag = 0 for linear and 1 for nonlinear derivative
329! (not needed for completely linear models)
330!
331!
332! Indices
333! x(i) res(j)
334! j: NL L=1
335! obj: NL
336!
337 k = 1
338 do j = 1, dimx
339 colsta(j) = k
340 do i = 1, nobs
341 rowno(k) = i
342 nlflag(k) = 1
343 k = k + 1
344 enddo
345 enddo
346 do i = 1, nobs
347 colsta(dimx+i) = k
348 rowno(k) = i
349 nlflag(k) = 0
350 value(k) = 1.d0
351 k = k + 1
352 rowno(k) = nobs+1
353 nlflag(k) = 1
354 k = k + 1
355 enddo
356 colsta(dimx+nobs+1) = k
357
358 lsq_readmatrix = 0 ! Return value means OK
359
360end Function lsq_readmatrix
361!
362!==========================================================================
363! Compute nonlinear terms and non-constant Jacobian elements
364!
365
366!> Compute nonlinear terms and non-constant Jacobian elements
367!!
368!! @include{doc} fdeval_params.dox
369Integer Function lsq_fdeval( x, g, jac, rowno, jcnm, mode, ignerr, errcnt, &
370 n, nz, thread, usrmem )
371#ifdef dec_directives_win32
372!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_FDEval
373#endif
374 use lsq_70
375 implicit none
376 integer, intent (in) :: n ! number of variables
377 integer, intent (in) :: rowno ! number of the row to be evaluated
378 integer, intent (in) :: nz ! number of nonzeros in this row
379 real*8, intent (in), dimension(n) :: x ! vector of current solution values
380 real*8, intent (in out) :: g ! constraint value
381 real*8, intent (in out), dimension(n) :: jac ! vector of derivatives for current constraint
382 integer, intent (in), dimension(nz) :: jcnm ! list of variables that appear nonlinearly
383 ! in this row. Ffor information only.
384 integer, intent (in) :: mode ! evaluation mode: 1 = function value
385 ! 2 = derivatives, 3 = both
386 integer, intent (in) :: ignerr ! if 1 then errors can be ignored as long
387 ! as errcnt is incremented
388 integer, intent (in out) :: errcnt ! error counter to be incremented in case
389 ! of function evaluation errors.
390 integer, intent (in) :: thread
391 real*8 usrmem(*) ! optional user memory
392
393 integer :: i,j
394 real*8 :: s
395!
396! Row Nobs+1: the objective function is nonlinear
397!
398 if ( rowno .eq. nobs+1 ) then
399!
400! Mode = 1 or 3. Function value: G = sum(i, res(i)**2)
401!
402 if ( mode .eq. 1 .or. mode .eq. 3 ) then
403 s = 0.d0
404 do i = 1, nobs
405 s = s + x(dimx+i)**2
406 enddo
407 g = s
408 endif
409!
410! Mode = 2 or 3: Derivative values:
411!
412 if ( mode .eq. 2 .or. mode .eq. 3 ) then
413 do i = 1, nobs
414 jac(dimx+i) = 2*x(dimx+i)
415 enddo
416 endif
417!
418! Row 1 to Nobs: The observation definitions:
419!
420 else
421!
422! Mode = 1 or 3: Function value - x-part only
423!
424 if ( mode .eq. 1 .or. mode .eq. 3 ) then
425 s = 0.d0
426 do j = 1, dimx
427 s = s + a(rowno,j)*x(j) + b(rowno,j)*x(j)**2
428 enddo
429 g = s
430 endif
431!
432! Mode = 2 or 3: Derivatives
433!
434 if ( mode .eq. 2 .or. mode .eq. 3 ) then
435 do j = 1, dimx
436 jac(j) = a(rowno,j) + 2.d0*b(rowno,j)*x(j)
437 enddo
438 endif
439 endif
440 lsq_fdeval = 0
441
442end Function lsq_fdeval
443
444
445!> Specify the structure of the Lagrangian of the Hessian
446!!
447!! @include{doc} 2DLagrStr_params.dox
448Integer Function lsq_2dlagrstr( HSRW, HSCL, NODRV, N, M, NHESS, UsrMem )
449#ifdef dec_directives_win32
450!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_2DLagrStr
451#endif
452 use lsq_70
453 Implicit None
454
455 Integer, Intent (IN) :: n, m, nhess
456 Integer, Intent (IN OUT) :: nodrv
457 Integer, Dimension(Nhess), Intent (Out) :: hsrw, hscl
458 real*8, Intent(IN OUT) :: usrmem(*)
459
460 integer :: i
461!
462! We assume that this is a Large Residual model, where the most
463! important 2nd derivatives are those that appear directly in the
464! objective.
465! We will therefore only define the 2nd derivatives corresponding to
466! the objective constraint where the 2nd derivatives is 2*I.
467! CONOPT will by default complain that some nonlinear variables do
468! not appear in the Hessian. We will therefore define one dummy element
469! on the diagonal of the Hessian for each of the nonlinear X-variables
470! and give it the numerical value 0.
471!
472! Because we use this simplitying assumption we must make sure that the
473! 2nd order information is not tested using option Lkdeb2 = 0
474!
475! Define structure of Hessian
476!
477 do i = 1, dimx+nobs
478 hsrw(i) = i
479 hscl(i) = i
480 enddo
481
482 lsq_2dlagrstr = 0
483
484End Function lsq_2dlagrstr
485
486
487!> Compute the Lagrangian of the Hessian
488!!
489!! @include{doc} 2DLagrVal_params.dox
490Integer Function lsq_2dlagrval( X, U, HSRW, HSCL, HSVL, NODRV, N, M, NHESS, UsrMem )
491#ifdef dec_directives_win32
492!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_2DLagrVal
493#endif
494 use lsq_70
495 Implicit None
496
497 Integer, Intent (IN) :: n, m, nhess
498 Integer, Intent (IN OUT) :: nodrv
499 real*8, Dimension(N), Intent (IN) :: x
500 real*8, Dimension(M), Intent (IN) :: u
501 Integer, Dimension(Nhess), Intent (In) :: hsrw, hscl
502 real*8, Dimension(NHess), Intent (Out) :: hsvl
503 real*8, Intent(IN OUT) :: usrmem(*)
504
505 integer :: i
506!
507! We assume that this is a Large Residual model, where the most
508! important 2nd derivatives are those that appear directly in the
509! objective.
510! We will therefore only define the 2nd derivatives corresponding to
511! the objective constraint where the 2nd derivatives is 2*I.
512! CONOPT will by default complain that some nonlinear variables do
513! not appear in the Hessian. We will therefore define one dummy element
514! on the diagonal of the Hessian for each of the nonlinear X-variables
515! and give it the numerical value 0.
516!
517! Normal Evaluation mode
518!
519 do i = 1, dimx
520 hsvl(i) = 0.d0
521 enddo
522 do i = 1, nobs
523 hsvl(dimx+i) = 2.d0*u(nobs+1)
524 enddo
525
526 lsq_2dlagrval = 0
527
528End Function lsq_2dlagrval
529
530
531!> Sets runtime options
532!!
533!! @include{doc} option_params.dox
534Integer Function lsq_option( ncall, rval, ival, lval, usrmem, name )
535#ifdef dec_directives_win32
536!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_Option
537#endif
538 integer ncall, ival, lval
539 character(Len=*) :: name
540 real*8 rval
541 real*8 usrmem(*) ! optional user memory
542
543 Select case (ncall)
544 case (1)
545 name = 'lkdeb2' ! Turn the debugger off
546 ival = 0
547 case default
548 name = ' '
549 end Select
550 lsq_option = 0
551
552end Function lsq_option
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
subroutine checkdual(case, minmax)
Definition comdecl.f90:394
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_option(cntvect, coi_option)
define callback routine for defining runtime options.
Definition conopt.f90:1346
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_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
integer function lsq_readmatrix(lower, curr, upper, vsta, type, rhs, esta, colsta, rowno, value, nlflag, n, m, nz, usrmem)
Define information about the model.
Definition leastl1.f90:166
integer function lsq_fdeval(x, g, jac, rowno, jcnm, mode, ignerr, errcnt, n, nz, thread, usrmem)
Compute nonlinear terms and non-constant Jacobian elements.
Definition leastl1.f90:287
integer function lsq_2dlagrstr(hsrw, hscl, nodrv, n, m, nhess, usrmem)
Specify the structure of the Lagrangian of the Hessian.
Definition leastsq2.f90:359
integer function lsq_2dlagrval(x, u, hsrw, hscl, hsvl, nodrv, n, m, nhess, usrmem)
Compute the Lagrangian of the Hessian.
Definition leastsq2.f90:396
integer function lsq_option(ncall, rval, ival, lval, usrmem, name)
Sets runtime options.
Definition leastsq2.f90:437
#define dimx
Definition leastsq.c:16
#define nobs
Definition leastsq.c:15
void defdata()
Defines the data for the problem.
Definition leastsq.c:35
float rndx()
Defines a pseudo random number between 0 and 1.
Definition leastsq.c:22
program leastsquare
Main program. A simple setup and call of CONOPT.
Definition leastsq.f90:68
real *8 obj
Definition comdecl.f90:16
integer solcalls
Definition comdecl.f90:15
integer sstat
Definition comdecl.f90:18
real *8, dimension(:), pointer udual
Definition comdecl.f90:24
integer, parameter minimize
Definition comdecl.f90:31
integer stacalls
Definition comdecl.f90:14
subroutine flog(msg, code)
Definition comdecl.f90:62
logical do_allocate
Definition comdecl.f90:27
real *8, dimension(:), pointer xprim
Definition comdecl.f90:23
integer mstat
Definition comdecl.f90:17
subroutine startup
Definition comdecl.f90:41