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)
169
170 call flog( "Successful Solve", 0 )
171!
172! Free solution memory
174 call finalize
175
176End Program leastsquare
177!
178! ============================================================================
179! Define information about the model:
180!
181
182!> Define information about the model
183!!
184!! @include{doc} readMatrix_params.dox
185Integer Function lsq_readmatrix( lower, curr, upper, vsta, type, rhs, esta, &
186 colsta, rowno, value, nlflag, n, m, nz, &
187 usrmem )
188#ifdef dec_directives_win32
189!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_ReadMatrix
190#endif
191 Use lsq_700
192 implicit none
193 integer, intent (in) :: n ! number of variables
194 integer, intent (in) :: m ! number of constraints
195 integer, intent (in) :: nz ! number of nonzeros
196 real*8, intent (in out), dimension(n) :: lower ! vector of lower bounds
197 real*8, intent (in out), dimension(n) :: curr ! vector of initial values
198 real*8, intent (in out), dimension(n) :: upper ! vector of upper bounds
199 integer, intent (in out), dimension(n) :: vsta ! vector of initial variable status
200 ! (not defined here)
201 integer, intent (out), dimension(m) :: type ! vector of equation types
202 integer, intent (in out), dimension(m) :: esta ! vector of initial equation status
203 ! (not defined here)
204 real*8, intent (in out), dimension(m) :: rhs ! vector of right hand sides
205 integer, intent (in out), dimension(n+1) :: colsta ! vector with start of column indices
206 integer, intent (out), dimension(nz) :: rowno ! vector of row numbers
207 integer, intent (in out), dimension(nz) :: nlflag ! vector of nonlinearity flags
208 real*8, intent (in out), dimension(nz) :: value ! vector of matrix values
209 real*8 usrmem(*) ! optional user memory
210
211 Integer :: i, j, k
212!
213! Information about Variables:
214! Default: Lower = -Inf, Curr = 0, and Upper = +inf.
215! Default: the status information in Vsta is not used.
216!
217 do i = 1, dimx
218 curr(i) = -0.8d0
219 enddo
220!
221! Information about Constraints:
222! Default: Rhs = 0
223! Default: the status information in Esta and the function
224! value in FV are not used.
225! Default: Type: There is no default.
226! 0 = Equality,
227! 1 = Greater than or equal,
228! 2 = Less than or equal,
229! 3 = Non binding.
230!
231! Constraints 1 to Nobs:
232! Rhs = Obs(i) and type Equality
233!
234 do i = 1, nobs
235 rhs(i) = obs(i)
236 type(i) = 0
237 enddo
238!
239! Constraint Nobs + 1 (Objective)
240! Rhs = 0 and type Non binding
241!
242 type(nobs+1) = 3
243!
244! Information about the Jacobian. CONOPT expects a columnwise
245! representation in Rowno, Value, Nlflag and Colsta.
246!
247! Colsta = Start of column indices (No Defaults):
248! Rowno = Row indices
249! Value = Value of derivative (by default only linear
250! derivatives are used)
251! Nlflag = 0 for linear and 1 for nonlinear derivative
252! (not needed for completely linear models)
253!
254!
255! Indices
256! x(i) res(j)
257! j: NL L=1
258! obj: NL
259! 3:
260!
261 k = 1
262 do j = 1, dimx
263 colsta(j) = k
264 do i = 1, nobs
265 rowno(k) = i
266 nlflag(k) = 1
267 k = k + 1
268 enddo
269 enddo
270 do i = 1, nobs
271 colsta(dimx+i) = k
272 rowno(k) = i
273 nlflag(k) = 0
274 value(k) = 1.d0
275 k = k + 1
276 rowno(k) = nobs+1
277 nlflag(k) = 1
278 k = k + 1
279 enddo
280 colsta(dimx+nobs+1) = k
281
282 lsq_readmatrix = 0 ! Return value means OK
283
284end Function lsq_readmatrix
285!
286!==========================================================================
287! Compute nonlinear terms and non-constant Jacobian elements
288!
289
290!> Compute nonlinear terms and non-constant Jacobian elements
291!!
292!! @include{doc} fdeval_params.dox
293Integer Recursive function lsq_fdeval( x, g, jac, rowno, jcnm, mode, ignerr, errcnt, &
294 n, nz, thread, usrmem )
295#ifdef dec_directives_win32
296!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_FDEval
297#endif
298 use lsq_700
299 implicit none
300 integer, intent (in) :: n ! number of variables
301 integer, intent (in) :: rowno ! number of the row to be evaluated
302 integer, intent (in) :: nz ! number of nonzeros in this row
303 real*8, intent (in), dimension(n) :: x ! vector of current solution values
304 real*8, intent (in out) :: g ! constraint value
305 real*8, intent (in out), dimension(n) :: jac ! vector of derivatives for current constraint
306 integer, intent (in), dimension(nz) :: jcnm ! list of variables that appear nonlinearly
307 ! in this row. Ffor information only.
308 integer, intent (in) :: mode ! evaluation mode: 1 = function value
309 ! 2 = derivatives, 3 = both
310 integer, intent (in) :: ignerr ! if 1 then errors can be ignored as long
311 ! as errcnt is incremented
312 integer, intent (in out) :: errcnt ! error counter to be incremented in case
313 ! of function evaluation errors.
314 integer, intent (in) :: thread
315 ! last call and 1 if the point is different.
316 real*8 usrmem(*) ! optional user memory
317
318 integer :: i,j
319 real*8 :: s
320!
321! Row 1: the objective function is nonlinear
322!
323 if ( rowno .eq. nobs+1 ) then
324!
325! Mode = 1 or 3. Function value: G = sum(i, res(i)**2)
326!
327 if ( mode .eq. 1 .or. mode .eq. 3 ) then
328 s = 0.d0
329 do i = 1, nobs
330 s = s + x(dimx+i)**2
331 enddo
332 g = s
333 endif
334!
335! Mode = 2 or 3: Derivative values:
336!
337 if ( mode .eq. 2 .or. mode .eq. 3 ) then
338 do i = 1, nobs
339 jac(dimx+i) = 2*x(dimx+i)
340 enddo
341 endif
342!
343! Row 1 to Nobs: The observation definitions:
344!
345 else
346!
347! Mode = 1 or 3: Function value - x-part only
348!
349 if ( mode .eq. 1 .or. mode .eq. 3 ) then
350 s = 0.d0
351 do j = 1, dimx
352 s = s + a(rowno,j)*x(j) + b(rowno,j)*x(j)**2
353 enddo
354 g = s
355 endif
356!
357! Mode = 2 or 3: Derivatives
358!
359 if ( mode .eq. 2 .or. mode .eq. 3 ) then
360 do j = 1, dimx
361 jac(j) = a(rowno,j) + 2.d0*b(rowno,j)*x(j)
362 enddo
363 endif
364 endif
365 lsq_fdeval = 0
366
367end Function lsq_fdeval
368
369
370!> Specify the structure of the Lagrangian of the Hessian
371!!
372!! @include{doc} 2DLagrStr_params.dox
373Integer Function lsq_2dlagrstr( HSRW, HSCL, NODRV, N, M, NHESS, UsrMem )
374#ifdef dec_directives_win32
375!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_2DLagrStr
376#endif
377 use lsq_700
378 Implicit None
379
380 Integer, Intent (IN) :: n, m, nhess
381 Integer, Intent (IN OUT) :: nodrv
382 Integer, Dimension(Nhess), Intent (Out) :: hsrw, hscl
383 real*8, Intent(IN OUT) :: usrmem(*)
384
385 integer :: i
386!
387! We assume that this is a Large Residual model, where the most
388! important 2nd derivatives are those that appear directly in the
389! objective.
390! We will therefore only define the 2nd derivatives corresponding to
391! the objective constraint where the 2nd derivatives is 2*I.
392! CONOPT will by default complain that some nonlinear variables do
393! not appear in the Hessian. We will therefore define one dummy element
394! on the diagonal of the Hessian for each of the nonlinear X-variables
395! and give it the numerical value 0.
396!
397!
398! Define structure of Hessian
399!
400 do i = 1, dimx+nobs
401 hsrw(i) = i
402 hscl(i) = i
403 enddo
404
405 lsq_2dlagrstr = 0
406
407End Function lsq_2dlagrstr
408
409
410!> Compute the Lagrangian of the Hessian
411!!
412!! @include{doc} 2DLagrVal_params.dox
413Integer Function lsq_2dlagrval( X, U, HSRW, HSCL, HSVL, NODRV, N, M, NHESS, UsrMem )
414#ifdef dec_directives_win32
415!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_2DLagrVal
416#endif
417 use lsq_700
418 Implicit None
419
420 Integer, Intent (IN) :: n, m, nhess
421 Integer, Intent (IN OUT) :: nodrv
422 real*8, Dimension(N), Intent (IN) :: x
423 real*8, Dimension(M), Intent (IN) :: u
424 Integer, Dimension(Nhess), Intent (In) :: hsrw, hscl
425 real*8, Dimension(NHess), Intent (Out) :: hsvl
426 real*8, Intent(IN OUT) :: usrmem(*)
427
428 integer :: i
429!
430! We assume that this is a Large Residual model, where the most
431! important 2nd derivatives are those that appear directly in the
432! objective.
433! We will therefore only define the 2nd derivatives corresponding to
434! the objective constraint where the 2nd derivatives is 2*I.
435! CONOPT will by default complain that some nonlinear variables do
436! not appear in the Hessian. We will therefore define one dummy element
437! on the diagonal of the Hessian for each of the nonlinear X-variables
438! and give it the numerical value 0.
439!
440! Normal Evaluation mode
441!
442 do i = 1, dimx
443 hsvl(i) = 0.d0
444 enddo
445 do i = 1, nobs
446 hsvl(dimx+i) = 2.d0*u(nobs+1)
447 enddo
448
449 lsq_2dlagrval = 0
450
451End Function lsq_2dlagrval
452
453
454!> Sets runtime options
455!!
456!! @include{doc} option_params.dox
457Integer Function lsq_option( ncall, rval, ival, lval, usrmem, name )
458#ifdef dec_directives_win32
459!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_Option
460#endif
461 integer ncall, ival, lval
462 character(Len=*) :: name
463 real*8 rval
464 real*8 usrmem(*) ! optional user memory
465
466 Select case (ncall)
467 case (1)
468 name = 'lkdeb2' ! Turn the debugger off
469 ival = 0
470 case default
471 name = ' '
472 end Select
473 lsq_option = 0
474
475end 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
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_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: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 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