CONOPT
Loading...
Searching...
No Matches
leastsq10.f90
Go to the documentation of this file.
1!> @file leastsq10.f90
2!! @ingroup FORT1THREAD_EXAMPLES
3!!
4!!
5!! This model is similar to leastsq5, but this time we define 2nd order
6!! directional information in the form of both an `2DDirIni` routine
7!! where the bulk of the work is done and a `2DDir` routine that just
8!! copies the information.
9!!
10!! We solve the following nonlinear least squares model:
11!!
12!! @verbatim
13!! min sum(i, res(i)**2 )
14!!
15!! sum(j, a(i,j)*x(j) + b(i,j)*x(j)**2 ) + res(i) = obs(i)
16!! @endverbatim
17!!
18!! where a, b, and obs are known data, and res and x are the
19!! variables of the model.
20!!
21!!
22!! For more information about the individual callbacks, please have a look at the source code.
23
24
25
26REAL FUNCTION rndx( )
27!
28! Defines a pseudo random number between 0 and 1
29!
30 IMPLICIT NONE
31
32 Integer, save :: seed = 12359
33
34 seed = mod(seed*1027+25,1048576)
35 rndx = float(seed)/float(1048576)
36
37END FUNCTION rndx
38
39subroutine defdata
40!
41! Define values for A, B, and Obs
42!
43 Use lsq_700
44 IMPLICIT NONE
45
46 Integer :: i, j
47 real*8, Parameter :: xtarg = -1.0
48 real*8, Parameter :: noise = 1.0
49 Real, External :: Rndx
50 real*8 :: o
51
52 do i = 1, nobs
53 o = 0.d0
54 do j = 1, dimx
55 a(i,j) = rndx()
56 b(i,j) = rndx()
57 o = o + a(i,j) * xtarg + b(i,j) * xtarg**2
58 enddo
59 obs(i) = o + noise * rndx()
60 enddo
61
62end subroutine defdata
63!
64! Main program.
65!
66!> Main program. A simple setup and call of CONOPT
67!!
69
70 Use proginfo
71 Use coidef
72 Use lsq_700
73 implicit None
74!
75! Declare the user callback routines as Integer, External:
76!
77 Integer, External :: lsq_readmatrix ! Mandatory Matrix definition routine defined below
78 Integer, External :: lsq_fdeval ! Function and Derivative evaluation routine
79 ! needed a nonlinear model.
80 Integer, External :: lsq_2ddirini ! Initialization routine for Directional 2nd derivatives
81 Integer, External :: lsq_2ddir ! Directional 2nd derivative routine
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#if defined(itl)
87!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_ReadMatrix
88!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_FDEval
89!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_2DDirIni
90!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_2DDir
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#endif
96!
97! Control vector
98!
99 INTEGER :: numcallback
100 INTEGER, Dimension(:), Pointer :: cntvect
101 INTEGER :: coi_error
102
103 call startup
104!
105! Define data
106!
107 Call defdata
108!
109! Create and initialize a Control Vector
110!
111 numcallback = coidef_size()
112 Allocate( cntvect(numcallback) )
113 coi_error = coidef_inifort( 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_optdir( cntvect, -1 ) ) ! Minimize
122 coi_error = max( coi_error, coidef_objcon( cntvect, nobs + 1 ) ) ! Objective is last constraint
123 coi_error = max( coi_error, coidef_optfile( cntvect, 'leastsq10.opt' ) )
124!
125! Tell CONOPT about the callback routines:
126!
127 coi_error = max( coi_error, coidef_readmatrix( cntvect, lsq_readmatrix ) )
128 coi_error = max( coi_error, coidef_fdeval( cntvect, lsq_fdeval ) )
129 coi_error = max( coi_error, coidef_2ddirini( cntvect, lsq_2ddirini ) )
130 coi_error = max( coi_error, coidef_2ddir( cntvect, lsq_2ddir ) )
131 coi_error = max( coi_error, coidef_status( cntvect, std_status ) )
132 coi_error = max( coi_error, coidef_solution( cntvect, std_solution ) )
133 coi_error = max( coi_error, coidef_message( cntvect, std_message ) )
134 coi_error = max( coi_error, coidef_errmsg( cntvect, std_errmsg ) )
135 coi_error = max( coi_error, coidef_debugfv( cntvect, 0 ) ) ! Debug Fdeval on or off
136
137#if defined(LICENSE_INT_1) && defined(LICENSE_INT_2) && defined(LICENSE_INT_3) && defined(LICENSE_TEXT)
138 coi_error = max( coi_error, coidef_license( cntvect, license_int_1, license_int_2, license_int_3, license_text) )
139#endif
140
141 If ( coi_error .ne. 0 ) THEN
142 write(*,*)
143 write(*,*) '**** Fatal Error while loading CONOPT Callback routines.'
144 write(*,*)
145 call flog( "Skipping Solve due to setup errors", 1 )
146 ENDIF
147!
148! Start CONOPT:
149!
150 coi_error = coi_solve( cntvect )
151
152 write(*,*)
153 write(*,*) 'End of Least Square example 10. Return code=',coi_error
154
155 If ( coi_error /= 0 ) then
156 call flog( "Errors encountered during solution", 1 )
157 elseif ( stacalls == 0 .or. solcalls == 0 ) then
158 call flog( "Status or Solution routine was not called", 1 )
159 elseif ( sstat /= 1 .or. mstat /= 2 ) then
160 call flog( "Solver and Model Status was not as expected (1,2)", 1 )
161 elseif ( abs( obj - 19.44434311d0 ) > 1.d-7 ) then
162 call flog( "Incorrect objective returned", 1 )
163 endif
164
165 if ( coi_free(cntvect) /= 0 ) call flog( "Error while freeing control vector",1)
166
167 call flog( "Successful Solve", 0 )
168
169End Program leastsquare
170!
171! ============================================================================
172! Define information about the model:
173!
174
175!> Define information about the model
176!!
177!! @include{doc} readMatrix_params.dox
178Integer Function lsq_readmatrix( lower, curr, upper, vsta, type, rhs, esta, &
179 colsta, rowno, value, nlflag, n, m, nz, &
180 usrmem )
181#if defined(itl)
182!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_ReadMatrix
183#endif
184 Use lsq_700
185 implicit none
186 integer, intent (in) :: n ! number of variables
187 integer, intent (in) :: m ! number of constraints
188 integer, intent (in) :: nz ! number of nonzeros
189 real*8, intent (in out), dimension(n) :: lower ! vector of lower bounds
190 real*8, intent (in out), dimension(n) :: curr ! vector of initial values
191 real*8, intent (in out), dimension(n) :: upper ! vector of upper bounds
192 integer, intent (in out), dimension(n) :: vsta ! vector of initial variable status
193 ! (not defined here)
194 integer, intent (out), dimension(m) :: type ! vector of equation types
195 integer, intent (in out), dimension(m) :: esta ! vector of initial equation status
196 ! (not defined here)
197 real*8, intent (in out), dimension(m) :: rhs ! vector of right hand sides
198 integer, intent (in out), dimension(n+1) :: colsta ! vector with start of column indices
199 integer, intent (out), dimension(nz) :: rowno ! vector of row numbers
200 integer, intent (in out), dimension(nz) :: nlflag ! vector of nonlinearity flags
201 real*8, intent (in out), dimension(nz) :: value ! vector of matrix values
202 real*8 usrmem(*) ! optional user memory
203
204 Integer :: i, j, k
205!
206! Information about Variables:
207! Default: Lower = -Inf, Curr = 0, and Upper = +inf.
208! Default: the status information in Vsta is not used.
209!
210 do i = 1, dimx
211 curr(i) = -0.8d0
212 enddo
213!
214! Information about Constraints:
215! Default: Rhs = 0
216! Default: the status information in Esta and the function
217! value in FV are not used.
218! Default: Type: There is no default.
219! 0 = Equality,
220! 1 = Greater than or equal,
221! 2 = Less than or equal,
222! 3 = Non binding.
223!
224! Constraints 1 to Nobs:
225! Rhs = Obs(i) and type Equality
226!
227 do i = 1, nobs
228 rhs(i) = obs(i)
229 type(i) = 0
230 enddo
231!
232! Constraint Nobs + 1 (Objective)
233! Rhs = 0 and type Non binding
234!
235 type(nobs+1) = 3
236!
237! Information about the Jacobian. We use the standard method with
238! Rowno, Value, Nlflag and Colsta and we do not use Colno.
239!
240! Colsta = Start of column indices (No Defaults):
241! Rowno = Row indices
242! Value = Value of derivative (by default only linear
243! derivatives are used)
244! Nlflag = 0 for linear and 1 for nonlinear derivative
245! (not needed for completely linear models)
246!
247!
248! Indices
249! x(i) res(j)
250! j: NL L=1
251! obj: NL
252! 3:
253!
254 k = 1
255 do j = 1, dimx
256 colsta(j) = k
257 do i = 1, nobs
258 rowno(k) = i
259 nlflag(k) = 1
260 k = k + 1
261 enddo
262 enddo
263 do i = 1, nobs
264 colsta(dimx+i) = k
265 rowno(k) = i
266 nlflag(k) = 0
267 value(k) = 1.d0
268 k = k + 1
269 rowno(k) = nobs+1
270 nlflag(k) = 1
271 k = k + 1
272 enddo
273 colsta(dimx+nobs+1) = k
274
275 lsq_readmatrix = 0 ! Return value means OK
276
277end Function lsq_readmatrix
278!
279!==========================================================================
280! Compute nonlinear terms and non-constant Jacobian elements
281!
282
283!> Compute nonlinear terms and non-constant Jacobian elements
284!!
285!! @include{doc} fdeval_params.dox
286Integer Function lsq_fdeval( x, g, jac, rowno, jcnm, mode, ignerr, errcnt, &
287 n, nz, thread, usrmem )
288#if defined(itl)
289!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_FDEval
290#endif
291 use lsq_700
292 implicit none
293 integer, intent (in) :: n ! number of variables
294 integer, intent (in) :: rowno ! number of the row to be evaluated
295 integer, intent (in) :: nz ! number of nonzeros in this row
296 real*8, intent (in), dimension(n) :: x ! vector of current solution values
297 real*8, intent (in out) :: g ! constraint value
298 real*8, intent (in out), dimension(n) :: jac ! vector of derivatives for current constraint
299 integer, intent (in), dimension(nz) :: jcnm ! list of variables that appear nonlinearly
300 ! in this row. Ffor information only.
301 integer, intent (in) :: mode ! evaluation mode: 1 = function value
302 ! 2 = derivatives, 3 = both
303 integer, intent (in) :: ignerr ! if 1 then errors can be ignored as long
304 ! as errcnt is incremented
305 integer, intent (in out) :: errcnt ! error counter to be incremented in case
306 ! of function evaluation errors.
307 integer, intent (in) :: thread
308 real*8 usrmem(*) ! optional user memory
309
310 integer :: i,j
311 real*8 :: s
312!
313! Row 1: the objective function is nonlinear
314!
315 if ( rowno .eq. nobs+1 ) then
316!
317! Mode = 1 or 3. Function value: G = sum(i, res(i)**2)
318!
319 if ( mode .eq. 1 .or. mode .eq. 3 ) then
320 s = 0.d0
321 do i = 1, nobs
322 s = s + x(dimx+i)**2
323 enddo
324 g = s
325 endif
326!
327! Mode = 2 or 3: Derivative values:
328!
329 if ( mode .eq. 2 .or. mode .eq. 3 ) then
330 do i = 1, nobs
331 jac(dimx+i) = 2*x(dimx+i)
332 enddo
333 endif
334!
335! Row 1 to Nobs: The observation definitions:
336!
337 else
338!
339! Mode = 1 or 3: Function value - x-part only
340!
341 if ( mode .eq. 1 .or. mode .eq. 3 ) then
342 s = 0.d0
343 do j = 1, dimx
344 s = s + a(rowno,j)*x(j) + b(rowno,j)*x(j)**2
345 enddo
346 g = s
347 endif
348!
349! Mode = 2 or 3: Derivatives
350!
351 if ( mode .eq. 2 .or. mode .eq. 3 ) then
352 do j = 1, dimx
353 jac(j) = a(rowno,j) + 2.d0*b(rowno,j)*x(j)
354 enddo
355 endif
356 endif
357 lsq_fdeval = 0
358
359end Function lsq_fdeval
360
361
362!> Computes the second derivative of a constraint in a direction
363!!
364!! @include{doc} 2DDir_params.dox
365Integer Function lsq_2ddir( X, DX, D2G, ROWNO, JCNM, NODRV, N, NJ, THREAD, USRMEM )
366#if defined(itl)
367!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_2DDir
368#endif
369 Use lsq_700
370 implicit none
371 INTEGER, Intent(IN) :: rowno, n, nj, thread
372 INTEGER, Intent(IN OUT) :: nodrv
373 INTEGER, Dimension(NJ), Intent(IN) :: jcnm
374 real*8, Dimension(N), Intent(IN) :: x
375 real*8, Dimension(N), Intent(IN) :: dx
376 real*8, Dimension(N), Intent(OUT) :: d2g
377 real*8, Intent(IN OUT) :: usrmem(*)
378 integer :: i, j
379
380 if ( rowno .eq. nobs+1 ) then
381!
382! Objective row: sum(i, res(i)**2 )
383!
384 do j = 1, dimx
385 d2g(j) = 0.d0
386 enddo
387 do i = 1, nobs
388 d2g(dimx+i) = 2.d0*dx(dimx+i)
389 enddo
390 else
391!
392! Constraint row with b(i,j)*x(j)**2 as only nonlinear term
393!
394 do j = 1, dimx
395 d2g(j) = c(rowno,j) ! Extract the value computed in 2DDirIni below
396 enddo
397 do i = 1, nobs
398 d2g(dimx+i) = 0.d0
399 enddo
400 endif
401
402 lsq_2ddir = 0
403
404End Function lsq_2ddir
405
406
407!> Computes the second derivative of a constraint in a direction
408!!
409!! @include{doc} 2DDir_params.dox
410
411!> Initialises the second derivative in a direction evaluation. Called each time the point or the direction changes
412!!
413!! @include{doc} 2DDirIni_params.dox
414Integer Function lsq_2ddirini( X, DX, RowList, ListSize, NumThread, NewPT, NODRV, NumVar, USRMEM )
415#if defined(itl)
416!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_2DDirIni
417#endif
418 Use lsq_700
419 implicit none
420 INTEGER, Intent(IN) :: listsize, numthread, newpt, numvar
421 Integer, Intent(IN OUT) :: nodrv
422 real*8, Dimension(NumVar), Intent(IN) :: x
423 real*8, Dimension(NumVar), Intent(IN) :: dx
424 Integer, Dimension(ListSize),Intent(IN) :: rowlist
425 real*8, Intent(IN OUT) :: usrmem(*)
426
427 integer :: i, j
428
429!
430! Constraint row i has b(i,j)*x(j)**2 as only nonlinear term and the
431! directional derivative is 2*b(i,j)*dx(j). The value is saved here and
432! used in the subsequent calls to Lsq_2DDir
433!
434 write(10,*) 'Inside 2DDirIni: NewPt=',newpt
435 Do i = 1, nobs
436 do j = 1, dimx
437 c(i,j) = 2.d0*b(i,j)*dx(j)
438 enddo
439 enddo
440
441 lsq_2ddirini = 0
442
443End Function lsq_2ddirini
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
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_2ddirini(cntvect, coi_2ddirini)
define callback routine for initializing the computation of second derivatives for a constraint in a ...
integer function coidef_2ddir(cntvect, coi_2ddir)
define callback routine for computing the second derivative for a constraint in a direction.
integer function coidef_optfile(cntvect, optfile)
define callback routine for defining an options file.
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_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_2ddirini(x, dx, rowlist, listsize, numthread, newpt, nodrv, numvar, usrmem)
Computes the second derivative of a constraint in a direction.
integer function lsq_2ddir(x, dx, d2g, rowno, jcnm, nodrv, n, nj, thread, usrmem)
Computes the second derivative of a constraint in a direction.
#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
#define nj
Definition mp_trans.c:46
real *8 obj
Definition comdecl.f90:10
integer solcalls
Definition comdecl.f90:9
integer sstat
Definition comdecl.f90:12
integer stacalls
Definition comdecl.f90:8
subroutine flog(msg, code)
Definition comdecl.f90:56
integer mstat
Definition comdecl.f90:11
subroutine startup
Definition comdecl.f90:35