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