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!! \f[
13!! \min \sum_i res_{i}^2 \\
14!!
15!! \sum_j ( a_{ij}x_j + b_{ij}x_j^2 ) + res_i = obs_i
16!! \f]
17!!
18!! where \f$a\f$, \f$b\f$, and \f$obs\f$ are known data, and \f$res\f$ and \f$x\f$ 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#if defined(_WIN32) && !defined(_WIN64)
25#define dec_directives_win32
26#endif
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!!
70Program leastsquare
71
73 Use conopt
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 :: lsq_2ddirini ! Initialization routine for Directional 2nd derivatives
83 Integer, External :: lsq_2ddir ! Directional 2nd derivative routine
84 Integer, External :: std_status ! Standard callback for displaying solution status
85 Integer, External :: std_solution ! Standard callback for displaying solution values
86 Integer, External :: std_message ! Standard callback for managing messages
87 Integer, External :: std_errmsg ! Standard callback for managing error messages
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 :: Lsq_2DDirIni
92!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_2DDir
93!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Status
94!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Solution
95!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Message
96!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_ErrMsg
97#endif
98!
99! Control vector
100!
101 INTEGER, Dimension(:), Pointer :: cntvect
102 INTEGER :: coi_error
103
104 call startup
105!
106! Define data
107!
108 Call defdata
109!
110! Create and initialize a Control Vector
111!
112 coi_error = coi_create( cntvect )
113!
114! Tell CONOPT about the size of the model by populating the Control Vector:
115!
116 coi_error = max( coi_error, coidef_numvar( cntvect, nobs + dimx ) )
117 coi_error = max( coi_error, coidef_numcon( cntvect, nobs + 1 ) )
118 coi_error = max( coi_error, coidef_numnz( cntvect, nobs * (dimx+2) ) )
119 coi_error = max( coi_error, coidef_numnlnz( cntvect, nobs * (dimx+1) ) )
120 coi_error = max( coi_error, coidef_optdir( cntvect, -1 ) ) ! Minimize
121 coi_error = max( coi_error, coidef_objcon( cntvect, nobs + 1 ) ) ! Objective is last constraint
122 coi_error = max( coi_error, coidef_optfile( cntvect, 'leastsq10.opt' ) )
123!
124! Tell CONOPT about the callback routines:
125!
126 coi_error = max( coi_error, coidef_readmatrix( cntvect, lsq_readmatrix ) )
127 coi_error = max( coi_error, coidef_fdeval( cntvect, lsq_fdeval ) )
128 coi_error = max( coi_error, coidef_2ddirini( cntvect, lsq_2ddirini ) )
129 coi_error = max( coi_error, coidef_2ddir( cntvect, lsq_2ddir ) )
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
136#if defined(CONOPT_LICENSE_INT_1) && defined(CONOPT_LICENSE_INT_2) && defined(CONOPT_LICENSE_INT_3) && defined(CONOPT_LICENSE_TEXT)
137 coi_error = max( coi_error, coidef_license( cntvect, conopt_license_int_1, conopt_license_int_2, conopt_license_int_3, conopt_license_text) )
138#endif
139
140 If ( coi_error .ne. 0 ) THEN
141 write(*,*)
142 write(*,*) '**** Fatal Error while loading CONOPT Callback routines.'
143 write(*,*)
144 call flog( "Skipping Solve due to setup errors", 1 )
145 ENDIF
146!
147! Start CONOPT:
148!
149 coi_error = coi_solve( cntvect )
150
151 write(*,*)
152 write(*,*) 'End of Least Square example 10. Return code=',coi_error
153
154 If ( coi_error /= 0 ) then
155 call flog( "Errors encountered during solution", 1 )
156 elseif ( stacalls == 0 .or. solcalls == 0 ) then
157 call flog( "Status or Solution routine was not called", 1 )
158 elseif ( sstat /= 1 .or. mstat /= 2 ) then
159 call flog( "Solver and Model Status was not as expected (1,2)", 1 )
160 elseif ( abs( obj - 19.44434311d0 ) > 1.d-7 ) then
161 call flog( "Incorrect objective returned", 1 )
162 endif
163
164 if ( coi_free(cntvect) /= 0 ) call flog( "Error while freeing control vector",1)
165
166 call flog( "Successful Solve", 0 )
167!
168! Free solution memory
169!
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 real*8 usrmem(*) ! optional user memory
312
313 integer :: i,j
314 real*8 :: s
315!
316! Row 1: the objective function is nonlinear
317!
318 if ( rowno .eq. nobs+1 ) then
319!
320! Mode = 1 or 3. Function value: G = sum(i, res(i)**2)
321!
322 if ( mode .eq. 1 .or. mode .eq. 3 ) then
323 s = 0.d0
324 do i = 1, nobs
325 s = s + x(dimx+i)**2
326 enddo
327 g = s
328 endif
329!
330! Mode = 2 or 3: Derivative values:
331!
332 if ( mode .eq. 2 .or. mode .eq. 3 ) then
333 do i = 1, nobs
334 jac(dimx+i) = 2*x(dimx+i)
335 enddo
336 endif
337!
338! Row 1 to Nobs: The observation definitions:
339!
340 else
341!
342! Mode = 1 or 3: Function value - x-part only
343!
344 if ( mode .eq. 1 .or. mode .eq. 3 ) then
345 s = 0.d0
346 do j = 1, dimx
347 s = s + a(rowno,j)*x(j) + b(rowno,j)*x(j)**2
348 enddo
349 g = s
350 endif
352! Mode = 2 or 3: Derivatives
353!
354 if ( mode .eq. 2 .or. mode .eq. 3 ) then
355 do j = 1, dimx
356 jac(j) = a(rowno,j) + 2.d0*b(rowno,j)*x(j)
357 enddo
358 endif
359 endif
360 lsq_fdeval = 0
361
362end Function lsq_fdeval
363
364
365!> Computes the second derivative of a constraint in a direction
366!!
367!! @include{doc} 2DDir_params.dox
368Integer Function lsq_2ddir( X, DX, D2G, ROWNO, JCNM, NODRV, N, NJ, THREAD, USRMEM )
369#ifdef dec_directives_win32
370!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_2DDir
371#endif
372 Use lsq_700
373 implicit none
374 INTEGER, Intent(IN) :: rowno, n, nj, thread
375 INTEGER, Intent(IN OUT) :: nodrv
376 INTEGER, Dimension(NJ), Intent(IN) :: jcnm
377 real*8, Dimension(N), Intent(IN) :: x
378 real*8, Dimension(N), Intent(IN) :: dx
379 real*8, Dimension(N), Intent(OUT) :: d2g
380 real*8, Intent(IN OUT) :: usrmem(*)
381 integer :: i, j
382
383 if ( rowno .eq. nobs+1 ) then
384!
385! Objective row: sum(i, res(i)**2 )
386!
387 do j = 1, dimx
388 d2g(j) = 0.d0
389 enddo
390 do i = 1, nobs
391 d2g(dimx+i) = 2.d0*dx(dimx+i)
392 enddo
393 else
394!
395! Constraint row with b(i,j)*x(j)**2 as only nonlinear term
396!
397 do j = 1, dimx
398 d2g(j) = c(rowno,j) ! Extract the value computed in 2DDirIni below
399 enddo
400 do i = 1, nobs
401 d2g(dimx+i) = 0.d0
402 enddo
403 endif
404
405 lsq_2ddir = 0
406
407End Function lsq_2ddir
408
409
410!> Computes the second derivative of a constraint in a direction
411!!
412!! @include{doc} 2DDir_params.dox
413
414!> Initialises the second derivative in a direction evaluation. Called each time the point or the direction changes
415!!
416!! @include{doc} 2DDirIni_params.dox
417Integer Function lsq_2ddirini( X, DX, RowList, ListSize, NumThread, NewPT, NODRV, NumVar, USRMEM )
418#ifdef dec_directives_win32
419!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_2DDirIni
420#endif
421 Use lsq_700
422 implicit none
423 INTEGER, Intent(IN) :: listsize, numthread, newpt, numvar
424 Integer, Intent(IN OUT) :: nodrv
425 real*8, Dimension(NumVar), Intent(IN) :: x
426 real*8, Dimension(NumVar), Intent(IN) :: dx
427 Integer, Dimension(ListSize),Intent(IN) :: rowlist
428 real*8, Intent(IN OUT) :: usrmem(*)
429
430 integer :: i, j
431
432!
433! Constraint row i has b(i,j)*x(j)**2 as only nonlinear term and the
434! directional derivative is 2*b(i,j)*dx(j). The value is saved here and
435! used in the subsequent calls to Lsq_2DDir
436!
437 write(10,*) 'Inside 2DDirIni: NewPt=',newpt
438 Do i = 1, nobs
439 do j = 1, dimx
440 c(i,j) = 2.d0*b(i,j)*dx(j)
441 enddo
442 enddo
443
444 lsq_2ddirini = 0
445
446End Function lsq_2ddirini
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_2ddirini(cntvect, coi_2ddirini)
define callback routine for initializing the computation of second derivatives for a constraint in a ...
Definition conopt.f90:1453
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_2ddir(cntvect, coi_2ddir)
define callback routine for computing the second derivative for a constraint in a direction.
Definition conopt.f90:1421
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_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_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: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
#define nj
Definition mp_trans.c:46
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