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