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!! \f[
14!! \min \sum_i res_{i}^2 \\
15!!
16!! \sum_j ( a_{ij}x_j + b_{ij}x_j^2 ) + res_i = obs_i
17!! \f]
18!!
19!! where \f$a\f$, \f$b\f$, and \f$obs\f$ are known data, and \f$res\f$ and \f$x\f$ are the
20!! variables of the model.
21!!
22!!
23!! For more information about the individual callbacks, please have a look at the source code.
24
25#if defined(_WIN32) && !defined(_WIN64)
26#define dec_directives_win32
27#endif
28
29module lsq7
30 Integer, parameter :: Nobs = 14
31 Integer, parameter :: DimX = 10
32 real*8, dimension(Nobs,DimX) :: a, b
33 real*8, dimension(Nobs) :: obs
34 Integer :: run ! Run = 1 standard,
35 ! Run = 2 with wide bounds on Res
36 ! Run = 3 standard model with obs(1) changed by +1.e-3
37end module lsq7
38
39
40REAL FUNCTION rndx( )
41!
42! Defines a pseudo random number between 0 and 1
43!
44 IMPLICIT NONE
45
46 Integer, save :: seed = 12359
47
48 seed = mod(seed*1027+25,1048576)
49 rndx = float(seed)/float(1048576)
50
51END FUNCTION rndx
52
53subroutine defdata
54!
55! Define values for A, B, and Obs
56!
57 Use lsq7
58 IMPLICIT NONE
59
60 Integer :: i, j
61 real*8, Parameter :: xtarg = -1.0
62 real*8, Parameter :: noise = 1.0
63 Real, External :: Rndx
64 real*8 :: o
65
66 do i = 1, nobs
67 o = 0.d0
68 do j = 1, dimx
69 a(i,j) = rndx()
70 b(i,j) = rndx()
71 o = o + a(i,j) * xtarg + b(i,j) * xtarg**2
72 enddo
73 obs(i) = o + noise * rndx()
74 enddo
75
76end subroutine defdata
77!
78! Main program.
79!
80!> Main program. A simple setup and call of CONOPT
81!!
82Program leastsquare
83
85 Use conopt
86 Use lsq7
87 implicit None
88!
89! Declare the user callback routines as Integer, External:
90!
91 Integer, External :: lsq_readmatrix ! Mandatory Matrix definition routine defined below
92 Integer, External :: lsq_fdeval ! Function and Derivative evaluation routine
93 ! needed a nonlinear model.
94 Integer, External :: std_status ! Standard callback for displaying solution status
95 Integer, External :: std_solution ! Standard callback for displaying solution values
96 Integer, External :: std_message ! Standard callback for managing messages
97 Integer, External :: std_errmsg ! Standard callback for managing error messages
98 Integer, External :: lsq_2dlagrstr ! Callback for second derivatives structure
99 Integer, External :: lsq_2dlagrval ! Callback for second derivatives values
100 Integer, External :: lsq_option ! Option defining callback routine
101#ifdef dec_directives_win32
102!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_ReadMatrix
103!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_FDEval
104!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Status
105!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Solution
106!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Message
107!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_ErrMsg
108!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_2DLagrStr
109!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_2DLagrVal
110!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_Option
111#endif
112!
113! Control vector
114!
115 INTEGER, Dimension(:), Pointer :: cntvect
116 INTEGER :: coi_error
117!
118 Integer :: i
119 Character(Len=132) :: text
120 real*8 :: obj_run1, obj_run3
121
122 call startup
123!
124! Define data
125!
126 Call defdata
127!
128! Create and initialize a Control Vector
129!
130 coi_error = coi_create( 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, 'leastsq7.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(CONOPT_LICENSE_INT_1) && defined(CONOPT_LICENSE_INT_2) && defined(CONOPT_LICENSE_INT_3) && defined(CONOPT_LICENSE_TEXT)
156 coi_error = max( coi_error, coidef_license( cntvect, conopt_license_int_1, conopt_license_int_2, conopt_license_int_3, conopt_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) ) > 1.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: Leastsq7', minimize )
191 endif
192 obj_run1 = obj
193 obj_run3 = obj_run1 + 1.0d-3*udual(1) ! expected new objective in run 3
194!
195! Repeat but with bounds on the residuals. They cannot be moved into
196! the post-triangle any more.
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) ) > 1.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: Leastsq7', minimize )
217 endif
218!
219! Repeat but with Obs(1) changed by 1.d-3. The objective should change
220! by 1.d-3*Udual(1).
221!
222 run = 3
223 obs(1) = obs(1) + 1.0d-3
224 write(10,*) 'Actual Objective for run 1=',obj_run1
225 write(10,*) 'Expected Objective for run 3=',obj_run3
226 coi_error = coi_solve( cntvect )
227 write(10,*) 'Actual Objective for run 3=',obj
228 write(10,*) 'Expected Objective for run 3=',obj_run3
229 write(10,*) 'Diffrens Objective for run 3=', abs(obj_run3-obj)
230
231 If ( coi_error /= 0 ) then
232 call flog( "Run 3: Errors encountered during solution", 1 )
233 elseif ( stacalls == 0 .or. solcalls == 0 ) then
234 call flog( "Run 3: Status or Solution routine was not called", 1 )
235 elseif ( .not. ( sstat == 1 .and. mstat == 2 ) ) then
236 call flog( "Run 3: Solver or Model status was not as expected (1,2)", 1 )
237 elseif ( abs( obj - obj_run3 ) > 1.d-2 ) then ! sloppy tolerance due to nonlinearity
238 call flog( "Run 3: Incorrect objective returned in run3 ", 1 )
239 else
240 Call checkdual( 'Run 3: Leastsq7', minimize )
241 endif
242
243 write(*,*)
244 write(*,*) 'End of Least Square example 7. Return code=',coi_error
245
246 if ( coi_free(cntvect) /= 0 ) call flog( "Error while freeing control vector",1)
247
248 call flog( "Successful Solve", 0 )
249!
250! Free solution memory
252 call finalize
253
254End Program leastsquare
255!
256! ============================================================================
257! Define information about the model:
258!
259
260!> Define information about the model
261!!
262!! @include{doc} readMatrix_params.dox
263Integer Function lsq_readmatrix( lower, curr, upper, vsta, type, rhs, esta, &
264 colsta, rowno, value, nlflag, n, m, nz, &
265 usrmem )
266#ifdef dec_directives_win32
267!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_ReadMatrix
268#endif
269 Use lsq7
270 implicit none
271 integer, intent (in) :: n ! number of variables
272 integer, intent (in) :: m ! number of constraints
273 integer, intent (in) :: nz ! number of nonzeros
274 real*8, intent (in out), dimension(n) :: lower ! vector of lower bounds
275 real*8, intent (in out), dimension(n) :: curr ! vector of initial values
276 real*8, intent (in out), dimension(n) :: upper ! vector of upper bounds
277 integer, intent (in out), dimension(n) :: vsta ! vector of initial variable status
278 ! (not defined here)
279 integer, intent (out), dimension(m) :: type ! vector of equation types
280 integer, intent (in out), dimension(m) :: esta ! vector of initial equation status
281 ! (not defined here)
282 real*8, intent (in out), dimension(m) :: rhs ! vector of right hand sides
283 integer, intent (in out), dimension(n+1) :: colsta ! vector with start of column indices
284 integer, intent (out), dimension(nz) :: rowno ! vector of row numbers
285 integer, intent (in out), dimension(nz) :: nlflag ! vector of nonlinearity flags
286 real*8, intent (in out), dimension(nz) :: value ! vector of matrix values
287 real*8 usrmem(*) ! optional user memory
288
289 Integer :: i, j, k
290!
291! Information about Variables:
292! Default: Lower = -Inf, Curr = 0, and Upper = +inf.
293! Default: the status information in Vsta is not used.
294!
295 do i = 1, dimx
296 curr(i) = -0.8d0
297 enddo
298 if ( run == 2 ) Then ! Add large bound on Res that will prevent the
299 do i = 1, nobs ! post triangle from being identified.
300 lower(dimx+i) = -1000.0d0
301 upper(dimx+i) = +1000.0d0
302 enddo
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. CONOPT expects a columnwise
329! representation in Rowno, Value, Nlflag and Colsta.
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#ifdef dec_directives_win32
379!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_FDEval
380#endif
381 use lsq7
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#ifdef dec_directives_win32
457!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_2DLagrStr
458#endif
459 use lsq7
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#ifdef dec_directives_win32
497!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_2DLagrVal
498#endif
499 use lsq7
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#ifdef dec_directives_win32
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: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 leastsq7.f90:34
integer run
Definition leastsq7.f90:36
real *8, dimension(nobs, dimx) b
Definition leastsq7.f90:34
real *8, dimension(nobs) obs
Definition leastsq7.f90:35
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