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)
243
244 call flog( "Successful Solve", 0 )
245!
246! Free solution memory
248 call finalize
249
250End Program leastsquare
251!
252! ============================================================================
253! Define information about the model:
254!
255
256!> Define information about the model
257!!
258!! @include{doc} readMatrix_params.dox
259Integer Function lsq_readmatrix( lower, curr, upper, vsta, type, rhs, esta, &
260 colsta, rowno, value, nlflag, n, m, nz, &
261 usrmem )
262#ifdef dec_directives_win32
263!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_ReadMatrix
264#endif
265 Use lsq_70
266 Use casedata_num
267 implicit none
268 integer, intent (in) :: n ! number of variables
269 integer, intent (in) :: m ! number of constraints
270 integer, intent (in) :: nz ! number of nonzeros
271 real*8, intent (in out), dimension(n) :: lower ! vector of lower bounds
272 real*8, intent (in out), dimension(n) :: curr ! vector of initial values
273 real*8, intent (in out), dimension(n) :: upper ! vector of upper bounds
274 integer, intent (in out), dimension(n) :: vsta ! vector of initial variable status
275 ! (not defined here)
276 integer, intent (out), dimension(m) :: type ! vector of equation types
277 integer, intent (in out), dimension(m) :: esta ! vector of initial equation status
278 ! (not defined here)
279 real*8, intent (in out), dimension(m) :: rhs ! vector of right hand sides
280 integer, intent (in out), dimension(n+1) :: colsta ! vector with start of column indices
281 integer, intent (out), dimension(nz) :: rowno ! vector of row numbers
282 integer, intent (in out), dimension(nz) :: nlflag ! vector of nonlinearity flags
283 real*8, intent (in out), dimension(nz) :: value ! vector of matrix values
284 real*8 usrmem(*) ! optional user memory
285
286 Integer :: i, j, k
287!
288! Information about Variables:
289! Default: Lower = -Inf, Curr = 0, and Upper = +inf.
290! Default: the status information in Vsta is not used.
291!
292 do i = 1, dimx
293 curr(i) = -0.8d0
294 enddo
295 if ( casenum == 2 ) Then ! Add large bound on Res that will prevent the
296 do i = 1, nobs ! post triangle from being identified.
297 lower(dimx+i) = -1000.0d0
298 upper(dimx+i) = +1000.0d0
299 enddo
300 endif
301!
302! Information about Constraints:
303! Default: Rhs = 0
304! Default: the status information in Esta and the function
305! value in FV are not used.
306! Default: Type: There is no default.
307! 0 = Equality,
308! 1 = Greater than or equal,
309! 2 = Less than or equal,
310! 3 = Non binding.
311!
312! Constraints 1 to Nobs:
313! Rhs = Obs(i) and type Equality
314!
315 do i = 1, nobs
316 rhs(i) = obs(i)
317 type(i) = 0
318 enddo
319!
320! Constraint Nobs + 1 (Objective)
321! Rhs = 0 and type Non binding
322!
323 type(nobs+1) = 3
324!
325! Information about the Jacobian. CONOPT expects a columnwise
326! representation in Rowno, Value, Nlflag and Colsta.
327!
328! Colsta = Start of column indices (No Defaults):
329! Rowno = Row indices
330! Value = Value of derivative (by default only linear
331! derivatives are used)
332! Nlflag = 0 for linear and 1 for nonlinear derivative
333! (not needed for completely linear models)
334!
335!
336! Indices
337! x(i) res(j)
338! j: NL L=1
339! obj: NL
340!
341 k = 1
342 do j = 1, dimx
343 colsta(j) = k
344 do i = 1, nobs
345 rowno(k) = i
346 nlflag(k) = 1
347 k = k + 1
348 enddo
349 enddo
350 do i = 1, nobs
351 colsta(dimx+i) = k
352 rowno(k) = i
353 nlflag(k) = 0
354 value(k) = 1.d0
355 k = k + 1
356 rowno(k) = nobs+1
357 nlflag(k) = 1
358 k = k + 1
359 enddo
360 colsta(dimx+nobs+1) = k
361
362 lsq_readmatrix = 0 ! Return value means OK
363
364end Function lsq_readmatrix
365!
366!==========================================================================
367! Compute nonlinear terms and non-constant Jacobian elements
368!
369
370!> Compute nonlinear terms and non-constant Jacobian elements
371!!
372!! @include{doc} fdeval_params.dox
373Integer Function lsq_fdeval( x, g, jac, rowno, jcnm, mode, ignerr, errcnt, &
374 n, nz, thread, usrmem )
375#ifdef dec_directives_win32
376!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_FDEval
377#endif
378 use lsq_70
379 implicit none
380 integer, intent (in) :: n ! number of variables
381 integer, intent (in) :: rowno ! number of the row to be evaluated
382 integer, intent (in) :: nz ! number of nonzeros in this row
383 real*8, intent (in), dimension(n) :: x ! vector of current solution values
384 real*8, intent (in out) :: g ! constraint value
385 real*8, intent (in out), dimension(n) :: jac ! vector of derivatives for current constraint
386 integer, intent (in), dimension(nz) :: jcnm ! list of variables that appear nonlinearly
387 ! in this row. Ffor information only.
388 integer, intent (in) :: mode ! evaluation mode: 1 = function value
389 ! 2 = derivatives, 3 = both
390 integer, intent (in) :: ignerr ! if 1 then errors can be ignored as long
391 ! as errcnt is incremented
392 integer, intent (in out) :: errcnt ! error counter to be incremented in case
393 ! of function evaluation errors.
394 integer, intent (in) :: thread
395 real*8 usrmem(*) ! optional user memory
396
397 integer :: i,j
398 real*8 :: s
399!
400! Row Nobs+1: the objective function is nonlinear
401!
402 if ( rowno .eq. nobs+1 ) then
403!
404! Mode = 1 or 3. Function value: G = sum(i, res(i)**2)
405!
406 if ( mode .eq. 1 .or. mode .eq. 3 ) then
407 s = 0.d0
408 do i = 1, nobs
409 s = s + x(dimx+i)**2
410 enddo
411 g = s
412 endif
413!
414! Mode = 2 or 3: Derivative values:
415!
416 if ( mode .eq. 2 .or. mode .eq. 3 ) then
417 do i = 1, nobs
418 jac(dimx+i) = 2*x(dimx+i)
419 enddo
420 endif
421!
422! Row 1 to Nobs: The observation definitions:
423!
424 else
425!
426! Mode = 1 or 3: Function value - x-part only
427!
428 if ( mode .eq. 1 .or. mode .eq. 3 ) then
429 s = 0.d0
430 do j = 1, dimx
431 s = s + a(rowno,j)*x(j) + b(rowno,j)*x(j)**2
432 enddo
433 g = s
434 endif
435!
436! Mode = 2 or 3: Derivatives
437!
438 if ( mode .eq. 2 .or. mode .eq. 3 ) then
439 do j = 1, dimx
440 jac(j) = a(rowno,j) + 2.d0*b(rowno,j)*x(j)
441 enddo
442 endif
443 endif
444 lsq_fdeval = 0
445
446end Function lsq_fdeval
447
448
449!> Specify the structure of the Lagrangian of the Hessian
450!!
451!! @include{doc} 2DLagrStr_params.dox
452Integer Function lsq_2dlagrstr( HSRW, HSCL, NODRV, N, M, NHESS, UsrMem )
453#ifdef dec_directives_win32
454!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_2DLagrStr
455#endif
456 use lsq_70
457 Implicit None
458
459 Integer, Intent (IN) :: n, m, nhess
460 Integer, Intent (IN OUT) :: nodrv
461 Integer, Dimension(Nhess), Intent (Out) :: hsrw, hscl
462 real*8, Intent(IN OUT) :: usrmem(*)
463
464 integer :: i
465!
466! We assume that this is a Large Residual model, where the most
467! important 2nd derivatives are those that appear directly in the
468! objective.
469! We will therefore only define the 2nd derivatives corresponding to
470! the objective constraint where the 2nd derivatives is 2*I.
471! CONOPT will by default complain that some nonlinear variables do
472! not appear in the Hessian. We will therefore define one dummy element
473! on the diagonal of the Hessian for each of the nonlinear X-variables
474! and give it the numerical value 0.
475!
476! Because we use this simplitying assumption we must make sure that the
477! 2nd order information is not tested using option Lkdeb2 = 0
478!
479! Define structure of Hessian
480!
481 do i = 1, dimx+nobs
482 hsrw(i) = i
483 hscl(i) = i
484 enddo
485
486 lsq_2dlagrstr = 0
487
488End Function lsq_2dlagrstr
489
490
491!> Compute the Lagrangian of the Hessian
492!!
493!! @include{doc} 2DLagrVal_params.dox
494Integer Function lsq_2dlagrval( X, U, HSRW, HSCL, HSVL, NODRV, N, M, NHESS, UsrMem )
495#ifdef dec_directives_win32
496!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_2DLagrVal
497#endif
498 use lsq_70
499 Implicit None
500
501 Integer, Intent (IN) :: n, m, nhess
502 Integer, Intent (IN OUT) :: nodrv
503 real*8, Dimension(N), Intent (IN) :: x
504 real*8, Dimension(M), Intent (IN) :: u
505 Integer, Dimension(Nhess), Intent (In) :: hsrw, hscl
506 real*8, Dimension(NHess), Intent (Out) :: hsvl
507 real*8, Intent(IN OUT) :: usrmem(*)
508
509 integer :: i
510!
511! We assume that this is a Large Residual model, where the most
512! important 2nd derivatives are those that appear directly in the
513! objective.
514! We will therefore only define the 2nd derivatives corresponding to
515! the objective constraint where the 2nd derivatives is 2*I.
516! CONOPT will by default complain that some nonlinear variables do
517! not appear in the Hessian. We will therefore define one dummy element
518! on the diagonal of the Hessian for each of the nonlinear X-variables
519! and give it the numerical value 0.
520!
521! Normal Evaluation mode
522!
523 do i = 1, dimx
524 hsvl(i) = 0.d0
525 enddo
526 do i = 1, nobs
527 hsvl(dimx+i) = 2.d0*u(nobs+1)
528 enddo
529
530 lsq_2dlagrval = 0
531
532End Function lsq_2dlagrval
533
534
535!> Sets runtime options
536!!
537!! @include{doc} option_params.dox
538Integer Function lsq_option( ncall, rval, ival, lval, usrmem, name )
539#ifdef dec_directives_win32
540!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_Option
541#endif
542 integer ncall, ival, lval
543 character(Len=*) :: name
544 real*8 rval
545 real*8 usrmem(*) ! optional user memory
546
547 Select case (ncall)
548 case (1)
549 name = 'lkdeb2' ! Turn the debugger off
550 ival = 0
551 case default
552 name = ' '
553 end Select
554 lsq_option = 0
555
556end Function lsq_option
integer function std_solution(xval, xmar, xbas, xsta, yval, ymar, ybas, ysta, n, m, usrmem)
Definition comdecl.f90:170
integer function std_status(modsta, solsta, iter, objval, usrmem)
Definition comdecl.f90:126
subroutine checkdual(case, minmax)
Definition comdecl.f90:432
integer function std_message(smsg, dmsg, nmsg, llen, usrmem, msgv)
Definition comdecl.f90:243
integer function std_errmsg(rowno, colno, posno, msglen, usrmem, msg)
Definition comdecl.f90:286
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:170
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:291
integer function lsq_2dlagrstr(hsrw, hscl, nodrv, n, m, nhess, usrmem)
Specify the structure of the Lagrangian of the Hessian.
Definition leastsq2.f90:363
integer function lsq_2dlagrval(x, u, hsrw, hscl, hsvl, nodrv, n, m, nhess, usrmem)
Compute the Lagrangian of the Hessian.
Definition leastsq2.f90:400
integer function lsq_option(ncall, rval, ival, lval, usrmem, name)
Sets runtime options.
Definition leastsq2.f90:441
#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
subroutine finalize
Definition comdecl.f90:79
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