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!! @verbatim
16!! min sum(i, res(i)**2 )
17!!
18!! sum(j, a(i,j)*x(j) + b(i,j)*x(j)**2 ) + res(i) = obs(i)
19!! @endverbatim
20!!
21!! where a, b, and obs are known data, and res and x 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
28REAL FUNCTION rndx( )
29!
30! Defines a pseudo random number between 0 and 1
31!
32 IMPLICIT NONE
33
34 Integer, save :: seed = 12359
35
36 seed = mod(seed*1027+25,1048576)
37 rndx = float(seed)/float(1048576)
38
39END FUNCTION rndx
40
41subroutine defdata
42!
43! Define values for A, B, and Obs
44!
45 Use lsq_700
46 IMPLICIT NONE
47
48 Integer :: i, j
49 real*8, Parameter :: xtarg = -1.0
50 real*8, Parameter :: noise = 1.0
51 Real, External :: Rndx
52 real*8 :: o
53
54 do i = 1, nobs
55 o = 0.d0
56 do j = 1, dimx
57 a(i,j) = rndx()
58 b(i,j) = rndx()
59 o = o + a(i,j) * xtarg + b(i,j) * xtarg**2
60 enddo
61 obs(i) = o + noise * rndx()
62 enddo
63
64end subroutine defdata
65!
66! Main program.
67!
68!> Main program. A simple setup and call of CONOPT
69!!
71
72 Use proginfo
73 Use coidef
74 Use lsq_700
75 implicit None
76!
77! Declare the user callback routines as Integer, External:
78!
79 Integer, External :: lsq_readmatrix ! Mandatory Matrix definition routine defined below
80 Integer, External :: lsq_fdeval ! Function and Derivative evaluation routine
81 ! needed a nonlinear model.
82 Integer, External :: std_status ! Standard callback for displaying solution status
83 Integer, External :: std_solution ! Standard callback for displaying solution values
84 Integer, External :: std_message ! Standard callback for managing messages
85 Integer, External :: std_errmsg ! Standard callback for managing error messages
86 Integer, External :: lsq_2dlagrstr ! Callback for second derivatives structure
87 Integer, External :: lsq_2dlagrval ! Callback for second derivatives values
88 Integer, External :: lsq_option ! Option defining callback routine
89#if defined(itl)
90!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_ReadMatrix
91!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_FDEval
92!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Status
93!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Solution
94!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Message
95!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_ErrMsg
96!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_2DLagrStr
97!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_2DLagrVal
98!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_Option
99#endif
100!
101! Control vector
102!
103 INTEGER :: numcallback
104 INTEGER, Dimension(:), Pointer :: cntvect
105 INTEGER :: coi_error
106
107 call startup
108!
109! Define data
110!
111 Call defdata
112!
113! Create and initialize a Control Vector
114!
115 numcallback = coidef_size()
116 Allocate( cntvect(numcallback) )
117 coi_error = coidef_inifort( 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(LICENSE_INT_1) && defined(LICENSE_INT_2) && defined(LICENSE_INT_3) && defined(LICENSE_TEXT)
144 coi_error = max( coi_error, coidef_license( cntvect, license_int_1, license_int_2, license_int_3, 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
181End Program leastsquare
182!
183! ============================================================================
184! Define information about the model:
185!
186
187!> Define information about the model
188!!
189!! @include{doc} readMatrix_params.dox
190Integer Function lsq_readmatrix( lower, curr, upper, vsta, type, rhs, esta, &
191 colsta, rowno, value, nlflag, n, m, nz, &
192 usrmem )
193#if defined(itl)
194!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_ReadMatrix
195#endif
196 Use lsq_700
197 implicit none
198 integer, intent (in) :: n ! number of variables
199 integer, intent (in) :: m ! number of constraints
200 integer, intent (in) :: nz ! number of nonzeros
201 real*8, intent (in out), dimension(n) :: lower ! vector of lower bounds
202 real*8, intent (in out), dimension(n) :: curr ! vector of initial values
203 real*8, intent (in out), dimension(n) :: upper ! vector of upper bounds
204 integer, intent (in out), dimension(n) :: vsta ! vector of initial variable status
205 ! (not defined here)
206 integer, intent (out), dimension(m) :: type ! vector of equation types
207 integer, intent (in out), dimension(m) :: esta ! vector of initial equation status
208 ! (not defined here)
209 real*8, intent (in out), dimension(m) :: rhs ! vector of right hand sides
210 integer, intent (in out), dimension(n+1) :: colsta ! vector with start of column indices
211 integer, intent (out), dimension(nz) :: rowno ! vector of row numbers
212 integer, intent (in out), dimension(nz) :: nlflag ! vector of nonlinearity flags
213 real*8, intent (in out), dimension(nz) :: value ! vector of matrix values
214 real*8 usrmem(*) ! optional user memory
215
216 Integer :: i, j, k
217!
218! Information about Variables:
219! Default: Lower = -Inf, Curr = 0, and Upper = +inf.
220! Default: the status information in Vsta is not used.
221!
222 do i = 1, dimx
223 curr(i) = -0.8d0
224 enddo
225!
226! Information about Constraints:
227! Default: Rhs = 0
228! Default: the status information in Esta and the function
229! value in FV are not used.
230! Default: Type: There is no default.
231! 0 = Equality,
232! 1 = Greater than or equal,
233! 2 = Less than or equal,
234! 3 = Non binding.
235!
236! Constraints 1 to Nobs:
237! Rhs = Obs(i) and type Equality
238!
239 do i = 1, nobs
240 rhs(i) = obs(i)
241 type(i) = 0
242 enddo
243!
244! Constraint Nobs + 1 (Objective)
245! Rhs = 0 and type Non binding
246!
247 type(nobs+1) = 3
248!
249! Information about the Jacobian. We use the standard method with
250! Rowno, Value, Nlflag and Colsta and we do not use Colno.
251!
252! Colsta = Start of column indices (No Defaults):
253! Rowno = Row indices
254! Value = Value of derivative (by default only linear
255! derivatives are used)
256! Nlflag = 0 for linear and 1 for nonlinear derivative
257! (not needed for completely linear models)
258!
259!
260! Indices
261! x(i) res(j)
262! j: NL L=1
263! obj: NL
264!
265 k = 1
266 do j = 1, dimx
267 colsta(j) = k
268 do i = 1, nobs
269 rowno(k) = i
270 nlflag(k) = 1
271 k = k + 1
272 enddo
273 enddo
274 do i = 1, nobs
275 colsta(dimx+i) = k
276 rowno(k) = i
277 nlflag(k) = 0
278 value(k) = 1.d0
279 k = k + 1
280 rowno(k) = nobs+1
281 nlflag(k) = 1
282 k = k + 1
283 enddo
284 colsta(dimx+nobs+1) = k
285
286 lsq_readmatrix = 0 ! Return value means OK
287
288end Function lsq_readmatrix
289!
290!==========================================================================
291! Compute nonlinear terms and non-constant Jacobian elements
292!
293
294!> Compute nonlinear terms and non-constant Jacobian elements
295!!
296!! @include{doc} fdeval_params.dox
297Integer Function lsq_fdeval( x, g, jac, rowno, jcnm, mode, ignerr, errcnt, &
298 n, nz, thread, usrmem )
299#if defined(itl)
300!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_FDEval
301#endif
302 use lsq_700
303 implicit none
304 integer, intent (in) :: n ! number of variables
305 integer, intent (in) :: rowno ! number of the row to be evaluated
306 integer, intent (in) :: nz ! number of nonzeros in this row
307 real*8, intent (in), dimension(n) :: x ! vector of current solution values
308 real*8, intent (in out) :: g ! constraint value
309 real*8, intent (in out), dimension(n) :: jac ! vector of derivatives for current constraint
310 integer, intent (in), dimension(nz) :: jcnm ! list of variables that appear nonlinearly
311 ! in this row. Ffor information only.
312 integer, intent (in) :: mode ! evaluation mode: 1 = function value
313 ! 2 = derivatives, 3 = both
314 integer, intent (in) :: ignerr ! if 1 then errors can be ignored as long
315 ! as errcnt is incremented
316 integer, intent (in out) :: errcnt ! error counter to be incremented in case
317 ! of function evaluation errors.
318 integer, intent (in) :: thread
319 real*8 usrmem(*) ! optional user memory
320
321 integer :: i,j
322 real*8 :: s
323!
324! Row Nobs+1: the objective function is nonlinear
325!
326 if ( rowno .eq. nobs+1 ) then
327!
328! Mode = 1 or 3. Function value: G = sum(i, res(i)**2)
329!
330 if ( mode .eq. 1 .or. mode .eq. 3 ) then
331 s = 0.d0
332 do i = 1, nobs
333 s = s + x(dimx+i)**2
334 enddo
335 g = s
336 endif
337!
338! Mode = 2 or 3: Derivative values:
339!
340 if ( mode .eq. 2 .or. mode .eq. 3 ) then
341 do i = 1, nobs
342 jac(dimx+i) = 2*x(dimx+i)
343 enddo
344 endif
345!
346! Row 1 to Nobs: The observation definitions:
347!
348 else
349!
350! Mode = 1 or 3: Function value - x-part only
351!
352 if ( mode .eq. 1 .or. mode .eq. 3 ) then
353 s = 0.d0
354 do j = 1, dimx
355 s = s + a(rowno,j)*x(j) + b(rowno,j)*x(j)**2
356 enddo
357 g = s
358 endif
359!
360! Mode = 2 or 3: Derivatives
361!
362 if ( mode .eq. 2 .or. mode .eq. 3 ) then
363 do j = 1, dimx
364 jac(j) = a(rowno,j) + 2.d0*b(rowno,j)*x(j)
365 enddo
366 endif
367 endif
368 lsq_fdeval = 0
369
370end Function lsq_fdeval
371
372
373!> Specify the structure of the Lagrangian of the Hessian
374!!
375!! @include{doc} 2DLagrStr_params.dox
376Integer Function lsq_2dlagrstr( HSRW, HSCL, NODRV, N, M, NHESS, UsrMem )
377#if defined(itl)
378!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_2DLagrStr
379#endif
380 use lsq_700
381 Implicit None
382
383 Integer, Intent (IN) :: n, m, nhess
384 Integer, Intent (IN OUT) :: nodrv
385 Integer, Dimension(Nhess), Intent (Out) :: hsrw, hscl
386 real*8, Intent(IN OUT) :: usrmem(*)
387
388 integer :: i
389!
390! We assume that this is a Large Residual model, where the most
391! important 2nd derivatives are those that appear directly in the
392! objective.
393! We will therefore only define the 2nd derivatives corresponding to
394! the objective constraint where the 2nd derivatives is 2*I.
395! CONOPT will by default complain that some nonlinear variables do
396! not appear in the Hessian. We will therefore define one dummy element
397! on the diagonal of the Hessian for each of the nonlinear X-variables
398! and give it the numerical value 0.
399!
400!
401! Define structure of Hessian
402!
403 do i = 1, dimx+nobs
404 hsrw(i) = i
405 hscl(i) = i
406 enddo
407
408 lsq_2dlagrstr = 0
409
410End Function lsq_2dlagrstr
411
412
413!> Compute the Lagrangian of the Hessian
414!!
415!! @include{doc} 2DLagrVal_params.dox
416Integer Function lsq_2dlagrval( X, U, HSRW, HSCL, HSVL, NODRV, N, M, NHESS, UsrMem )
417#if defined(itl)
418!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_2DLagrVal
419#endif
420 use lsq_700
421 Implicit None
422
423 Integer, Intent (IN) :: n, m, nhess
424 Integer, Intent (IN OUT) :: nodrv
425 real*8, Dimension(N), Intent (IN) :: x
426 real*8, Dimension(M), Intent (IN) :: u
427 Integer, Dimension(Nhess), Intent (In) :: hsrw, hscl
428 real*8, Dimension(NHess), Intent (Out) :: hsvl
429 real*8, Intent(IN OUT) :: usrmem(*)
430
431 integer :: i
432!
433! We assume that this is a Large Residual model, where the most
434! important 2nd derivatives are those that appear directly in the
435! objective.
436! We will therefore only define the 2nd derivatives corresponding to
437! the objective constraint where the 2nd derivatives is 2*I.
438! CONOPT will by default complain that some nonlinear variables do
439! not appear in the Hessian. We will therefore define one dummy element
440! on the diagonal of the Hessian for each of the nonlinear X-variables
441! and give it the numerical value 0.
442!
443! Normal Evaluation mode
444!
445 do i = 1, dimx
446 hsvl(i) = 0.d0
447 enddo
448 do i = 1, nobs
449 hsvl(dimx+i) = 2.d0*u(nobs+1)
450 enddo
451
452 lsq_2dlagrval = 0
453
454End Function lsq_2dlagrval
455
456
457!> Sets runtime options
458!!
459!! @include{doc} option_params.dox
460Integer Function lsq_option( ncall, rval, ival, lval, usrmem, name )
461#if defined(itl)
462!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_Option
463#endif
464 integer ncall, ival, lval
465 character(Len=*) :: name
466 real*8 rval
467 real*8 usrmem(*) ! optional user memory
468
469 Select case (ncall)
470 case (1)
471 name = 'lkdeb2' ! Turn the debugger off
472 ival = 0
473 case default
474 name = ' '
475 end Select
476 lsq_option = 0
477
478end 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_debugfv(cntvect, debugfv)
turn Debugging of FDEval on and off.
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 obj
Definition comdecl.f90:10
integer solcalls
Definition comdecl.f90:9
integer sstat
Definition comdecl.f90:12
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
integer mstat
Definition comdecl.f90:11
subroutine startup
Definition comdecl.f90:35