CONOPT
Loading...
Searching...
No Matches
leastsq2.f90
Go to the documentation of this file.
1!> @file leastsq2.f90
2!! @ingroup FORT1THREAD_EXAMPLES
3!!
4!!
5!! This model is similar to leastsq. The key difference is that
6!! we supply a callback routine that can compute 2nd derivatives of
7!! the model. However, we only include part of the 2nd derivatives
8!! corresponding to the direct objective terms, res(i)**2.
9!! The terms from b(i,j)*x(j)**2 are ignored in the 2nd derivatives.
10!! CONOPT will not notice the incorrect derivatives but it may converge
11!! more slowly.
12!!
13!! We solve the following nonlinear least squares model:
14!!
15!! \f[
16!! \min \sum_i res_{i}^2 \\
17!!
18!! \sum_j ( a_{ij}x_j + b_{ij}x_j^2 ) + res_i = obs_i
19!! \f]
20!!
21!! where \f$a\f$, \f$b\f$, and \f$obs\f$ are known data, and \f$res\f$ and \f$x\f$ are the
22!! variables of the model.
23!!
24!!
25!! For more information about the individual callbacks, please have a look at the source code.
26
27#if defined(_WIN32) && !defined(_WIN64)
28#define dec_directives_win32
29#endif
30
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_700
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!!
73Program leastsquare
74
76 Use conopt
77 Use lsq_700
78 implicit None
79!
80! Declare the user callback routines as Integer, External:
81!
82 Integer, External :: lsq_readmatrix ! Mandatory Matrix definition routine defined below
83 Integer, External :: lsq_fdeval ! Function and Derivative evaluation routine
84 ! needed a nonlinear model.
85 Integer, External :: std_status ! Standard callback for displaying solution status
86 Integer, External :: std_solution ! Standard callback for displaying solution values
87 Integer, External :: std_message ! Standard callback for managing messages
88 Integer, External :: std_errmsg ! Standard callback for managing error messages
89 Integer, External :: lsq_2dlagrstr ! Callback for second derivatives structure
90 Integer, External :: lsq_2dlagrval ! Callback for second derivatives values
91 Integer, External :: lsq_option ! Option defining callback routine
92#ifdef dec_directives_win32
93!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_ReadMatrix
94!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_FDEval
95!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Status
96!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Solution
97!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Message
98!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_ErrMsg
99!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_2DLagrStr
100!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_2DLagrVal
101!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_Option
102#endif
103!
104! Control vector
105!
106 INTEGER, Dimension(:), Pointer :: cntvect
107 INTEGER :: coi_error
108
109 call startup
110!
111! Define data
112!
113 Call defdata
114!
115! Create and initialize a Control Vector
116!
117 coi_error = coi_create( cntvect )
118!
119! Tell CONOPT about the size of the model by populating the Control Vector:
120!
121 coi_error = max( coi_error, coidef_numvar( cntvect, nobs + dimx ) )
122 coi_error = max( coi_error, coidef_numcon( cntvect, nobs + 1 ) )
123 coi_error = max( coi_error, coidef_numnz( cntvect, nobs * (dimx+2) ) )
124 coi_error = max( coi_error, coidef_numnlnz( cntvect, nobs * (dimx+1) ) )
125 coi_error = max( coi_error, coidef_numhess( cntvect, nobs + dimx ) )
126 coi_error = max( coi_error, coidef_optdir( cntvect, -1 ) ) ! Minimize
127 coi_error = max( coi_error, coidef_objcon( cntvect, nobs + 1 ) ) ! Objective is last constraint
128 coi_error = max( coi_error, coidef_debugfv( cntvect, 0 ) ) ! Debug Fdeval on or off
129 coi_error = max( coi_error, coidef_optfile( cntvect, 'leastsq2.opt' ) )
130!
131! Tell CONOPT about the callback routines:
132!
133 coi_error = max( coi_error, coidef_readmatrix( cntvect, lsq_readmatrix ) )
134 coi_error = max( coi_error, coidef_fdeval( cntvect, lsq_fdeval ) )
135 coi_error = max( coi_error, coidef_status( cntvect, std_status ) )
136 coi_error = max( coi_error, coidef_solution( cntvect, std_solution ) )
137 coi_error = max( coi_error, coidef_message( cntvect, std_message ) )
138 coi_error = max( coi_error, coidef_errmsg( cntvect, std_errmsg ) )
139 coi_error = max( coi_error, coidef_2dlagrstr( cntvect, lsq_2dlagrstr ) )
140 coi_error = max( coi_error, coidef_2dlagrval( cntvect, lsq_2dlagrval ) )
141 coi_error = max( coi_error, coidef_option( cntvect, lsq_option ) )
142
143#if defined(CONOPT_LICENSE_INT_1) && defined(CONOPT_LICENSE_INT_2) && defined(CONOPT_LICENSE_INT_3) && defined(CONOPT_LICENSE_TEXT)
144 coi_error = max( coi_error, coidef_license( cntvect, conopt_license_int_1, conopt_license_int_2, conopt_license_int_3, conopt_license_text) )
145#endif
146
147 If ( coi_error .ne. 0 ) THEN
148 write(*,*)
149 write(*,*) '**** Fatal Error while loading CONOPT Callback routines.'
150 write(*,*)
151 call flog( "Skipping Solve due to setup errors", 1 )
152 ENDIF
153!
154! Save the solution so we can check the duals:
155!
156 do_allocate = .true.
157!
158! Start CONOPT:
159!
160 coi_error = coi_solve( cntvect )
161
162 write(*,*)
163 write(*,*) 'End of Least Square example 2. Return code=',coi_error
164
165 If ( coi_error /= 0 ) then
166 call flog( "Errors encountered during solution", 1 )
167 elseif ( stacalls == 0 .or. solcalls == 0 ) then
168 call flog( "Status or Solution routine was not called", 1 )
169 elseif ( sstat /= 1 .or. mstat /= 2 ) then
170 call flog( "Solver or Model status was not as expected (1,2)", 1 )
171 elseif ( abs( obj - 19.44434311d0 ) > 1.d-7 ) then
172 call flog( "Incorrect objective returned", 1 )
173 Else
174 Call checkdual( 'LeastSq2', minimize )
175 endif
176
177 if ( coi_free(cntvect) /= 0 ) call flog( "Error while freeing control vector",1)
178
179 call flog( "Successful Solve", 0 )
180!
181! Free solution memory
183 call finalize
184
185End Program leastsquare
186!
187! ============================================================================
188! Define information about the model:
189!
190
191!> Define information about the model
192!!
193!! @include{doc} readMatrix_params.dox
194Integer Function lsq_readmatrix( lower, curr, upper, vsta, type, rhs, esta, &
195 colsta, rowno, value, nlflag, n, m, nz, &
196 usrmem )
197#ifdef dec_directives_win32
198!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_ReadMatrix
199#endif
200 Use lsq_700
201 implicit none
202 integer, intent (in) :: n ! number of variables
203 integer, intent (in) :: m ! number of constraints
204 integer, intent (in) :: nz ! number of nonzeros
205 real*8, intent (in out), dimension(n) :: lower ! vector of lower bounds
206 real*8, intent (in out), dimension(n) :: curr ! vector of initial values
207 real*8, intent (in out), dimension(n) :: upper ! vector of upper bounds
208 integer, intent (in out), dimension(n) :: vsta ! vector of initial variable status
209 ! (not defined here)
210 integer, intent (out), dimension(m) :: type ! vector of equation types
211 integer, intent (in out), dimension(m) :: esta ! vector of initial equation status
212 ! (not defined here)
213 real*8, intent (in out), dimension(m) :: rhs ! vector of right hand sides
214 integer, intent (in out), dimension(n+1) :: colsta ! vector with start of column indices
215 integer, intent (out), dimension(nz) :: rowno ! vector of row numbers
216 integer, intent (in out), dimension(nz) :: nlflag ! vector of nonlinearity flags
217 real*8, intent (in out), dimension(nz) :: value ! vector of matrix values
218 real*8 usrmem(*) ! optional user memory
219
220 Integer :: i, j, k
221!
222! Information about Variables:
223! Default: Lower = -Inf, Curr = 0, and Upper = +inf.
224! Default: the status information in Vsta is not used.
225!
226 do i = 1, dimx
227 curr(i) = -0.8d0
228 enddo
229!
230! Information about Constraints:
231! Default: Rhs = 0
232! Default: the status information in Esta and the function
233! value in FV are not used.
234! Default: Type: There is no default.
235! 0 = Equality,
236! 1 = Greater than or equal,
237! 2 = Less than or equal,
238! 3 = Non binding.
239!
240! Constraints 1 to Nobs:
241! Rhs = Obs(i) and type Equality
242!
243 do i = 1, nobs
244 rhs(i) = obs(i)
245 type(i) = 0
246 enddo
247!
248! Constraint Nobs + 1 (Objective)
249! Rhs = 0 and type Non binding
250!
251 type(nobs+1) = 3
252!
253! Information about the Jacobian. CONOPT expects a columnwise
254! representation in Rowno, Value, Nlflag and Colsta.
255!
256! Colsta = Start of column indices (No Defaults):
257! Rowno = Row indices
258! Value = Value of derivative (by default only linear
259! derivatives are used)
260! Nlflag = 0 for linear and 1 for nonlinear derivative
261! (not needed for completely linear models)
262!
263!
264! Indices
265! x(i) res(j)
266! j: NL L=1
267! obj: NL
268!
269 k = 1
270 do j = 1, dimx
271 colsta(j) = k
272 do i = 1, nobs
273 rowno(k) = i
274 nlflag(k) = 1
275 k = k + 1
276 enddo
277 enddo
278 do i = 1, nobs
279 colsta(dimx+i) = k
280 rowno(k) = i
281 nlflag(k) = 0
282 value(k) = 1.d0
283 k = k + 1
284 rowno(k) = nobs+1
285 nlflag(k) = 1
286 k = k + 1
287 enddo
288 colsta(dimx+nobs+1) = k
289
290 lsq_readmatrix = 0 ! Return value means OK
291
292end Function lsq_readmatrix
293!
294!==========================================================================
295! Compute nonlinear terms and non-constant Jacobian elements
296!
297
298!> Compute nonlinear terms and non-constant Jacobian elements
299!!
300!! @include{doc} fdeval_params.dox
301Integer Recursive function lsq_fdeval( x, g, jac, rowno, jcnm, mode, ignerr, errcnt, &
302 n, nz, thread, usrmem )
303#ifdef dec_directives_win32
304!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_FDEval
305#endif
306 use lsq_700
307 implicit none
308 integer, intent (in) :: n ! number of variables
309 integer, intent (in) :: rowno ! number of the row to be evaluated
310 integer, intent (in) :: nz ! number of nonzeros in this row
311 real*8, intent (in), dimension(n) :: x ! vector of current solution values
312 real*8, intent (in out) :: g ! constraint value
313 real*8, intent (in out), dimension(n) :: jac ! vector of derivatives for current constraint
314 integer, intent (in), dimension(nz) :: jcnm ! list of variables that appear nonlinearly
315 ! in this row. Ffor information only.
316 integer, intent (in) :: mode ! evaluation mode: 1 = function value
317 ! 2 = derivatives, 3 = both
318 integer, intent (in) :: ignerr ! if 1 then errors can be ignored as long
319 ! as errcnt is incremented
320 integer, intent (in out) :: errcnt ! error counter to be incremented in case
321 ! of function evaluation errors.
322 integer, intent (in) :: thread
323 real*8 usrmem(*) ! optional user memory
324
325 integer :: i,j
326 real*8 :: s
327!
328! Row Nobs+1: the objective function is nonlinear
329!
330 if ( rowno .eq. nobs+1 ) then
331!
332! Mode = 1 or 3. Function value: G = sum(i, res(i)**2)
333!
334 if ( mode .eq. 1 .or. mode .eq. 3 ) then
335 s = 0.d0
336 do i = 1, nobs
337 s = s + x(dimx+i)**2
338 enddo
339 g = s
340 endif
341!
342! Mode = 2 or 3: Derivative values:
343!
344 if ( mode .eq. 2 .or. mode .eq. 3 ) then
345 do i = 1, nobs
346 jac(dimx+i) = 2*x(dimx+i)
347 enddo
348 endif
349!
350! Row 1 to Nobs: The observation definitions:
351!
352 else
353!
354! Mode = 1 or 3: Function value - x-part only
355!
356 if ( mode .eq. 1 .or. mode .eq. 3 ) then
357 s = 0.d0
358 do j = 1, dimx
359 s = s + a(rowno,j)*x(j) + b(rowno,j)*x(j)**2
360 enddo
361 g = s
362 endif
363!
364! Mode = 2 or 3: Derivatives
365!
366 if ( mode .eq. 2 .or. mode .eq. 3 ) then
367 do j = 1, dimx
368 jac(j) = a(rowno,j) + 2.d0*b(rowno,j)*x(j)
369 enddo
370 endif
371 endif
372 lsq_fdeval = 0
373
374end Function lsq_fdeval
375
376
377!> Specify the structure of the Lagrangian of the Hessian
378!!
379!! @include{doc} 2DLagrStr_params.dox
380Integer Function lsq_2dlagrstr( HSRW, HSCL, NODRV, N, M, NHESS, UsrMem )
381#ifdef dec_directives_win32
382!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_2DLagrStr
383#endif
384 use lsq_700
385 Implicit None
386
387 Integer, Intent (IN) :: n, m, nhess
388 Integer, Intent (IN OUT) :: nodrv
389 Integer, Dimension(Nhess), Intent (Out) :: hsrw, hscl
390 real*8, Intent(IN OUT) :: usrmem(*)
391
392 integer :: i
393!
394! We assume that this is a Large Residual model, where the most
395! important 2nd derivatives are those that appear directly in the
396! objective.
397! We will therefore only define the 2nd derivatives corresponding to
398! the objective constraint where the 2nd derivatives is 2*I.
399! CONOPT will by default complain that some nonlinear variables do
400! not appear in the Hessian. We will therefore define one dummy element
401! on the diagonal of the Hessian for each of the nonlinear X-variables
402! and give it the numerical value 0.
403!
404!
405! Define structure of Hessian
406!
407 do i = 1, dimx+nobs
408 hsrw(i) = i
409 hscl(i) = i
410 enddo
411
412 lsq_2dlagrstr = 0
413
414End Function lsq_2dlagrstr
415
416
417!> Compute the Lagrangian of the Hessian
418!!
419!! @include{doc} 2DLagrVal_params.dox
420Integer Function lsq_2dlagrval( X, U, HSRW, HSCL, HSVL, NODRV, N, M, NHESS, UsrMem )
421#ifdef dec_directives_win32
422!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_2DLagrVal
423#endif
424 use lsq_700
425 Implicit None
426
427 Integer, Intent (IN) :: n, m, nhess
428 Integer, Intent (IN OUT) :: nodrv
429 real*8, Dimension(N), Intent (IN) :: x
430 real*8, Dimension(M), Intent (IN) :: u
431 Integer, Dimension(Nhess), Intent (In) :: hsrw, hscl
432 real*8, Dimension(NHess), Intent (Out) :: hsvl
433 real*8, Intent(IN OUT) :: usrmem(*)
434
435 integer :: i
436!
437! We assume that this is a Large Residual model, where the most
438! important 2nd derivatives are those that appear directly in the
439! objective.
440! We will therefore only define the 2nd derivatives corresponding to
441! the objective constraint where the 2nd derivatives is 2*I.
442! CONOPT will by default complain that some nonlinear variables do
443! not appear in the Hessian. We will therefore define one dummy element
444! on the diagonal of the Hessian for each of the nonlinear X-variables
445! and give it the numerical value 0.
446!
447! Normal Evaluation mode
448!
449 do i = 1, dimx
450 hsvl(i) = 0.d0
451 enddo
452 do i = 1, nobs
453 hsvl(dimx+i) = 2.d0*u(nobs+1)
454 enddo
455
456 lsq_2dlagrval = 0
457
458End Function lsq_2dlagrval
459
460
461!> Sets runtime options
462!!
463!! @include{doc} option_params.dox
464Integer Function lsq_option( ncall, rval, ival, lval, usrmem, name )
465#ifdef dec_directives_win32
466!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_Option
467#endif
468 integer ncall, ival, lval
469 character(Len=*) :: name
470 real*8 rval
471 real*8 usrmem(*) ! optional user memory
472
473 Select case (ncall)
474 case (1)
475 name = 'lkdeb2' ! Turn the debugger off
476 ival = 0
477 case default
478 name = ' '
479 end Select
480 lsq_option = 0
481
482end 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_debugfv(cntvect, debugfv)
turn Debugging of FDEval on and off.
Definition conopt.f90:387
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 obj
Definition comdecl.f90:16
integer solcalls
Definition comdecl.f90:15
integer sstat
Definition comdecl.f90:18
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
integer mstat
Definition comdecl.f90:17
subroutine startup
Definition comdecl.f90:41