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