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