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