CONOPT
Loading...
Searching...
No Matches
leastsq.f90
Go to the documentation of this file.
1!> @file leastsq.f90
2!! @ingroup FORT1THREAD_EXAMPLES
3!!
4!! Solves a nonlinear least squares model
5!!
6!! We solve the following nonlinear least squares model:
7!!
8!! \f[
9!! \min \sum_i res_{i}^2 \\
10!!
11!! \sum_j ( a_{ij}x_j + b_{ij}x_j^2 ) + res_i = obs_i
12!! \f]
13!!
14!! where \f$a\f$, \f$b\f$, and \f$obs\f$ are known data, and \f$res\f$ and \f$x\f$ are the
15!! variables of the model.
16!!
17!!
18!! For more information about the individual callbacks, please have a look at the source code.
19
20#if defined(_WIN32) && !defined(_WIN64)
21#define dec_directives_win32
22#endif
23
24REAL FUNCTION rndx( )
25!
26! Defines a pseudo random number between 0 and 1
27!
28 IMPLICIT NONE
29
30 Integer, save :: seed = 12359
31
32 seed = mod(seed*1027+25,1048576)
33 rndx = float(seed)/float(1048576)
34
35END FUNCTION rndx
36
37subroutine defdata
38!
39! Define values for A, B, and Obs
40!
41 Use lsq_700
42 IMPLICIT NONE
43
44 Integer :: i, j
45 real*8, Parameter :: xtarg = -1.0
46 real*8, Parameter :: noise = 1.0
47 Real, External :: Rndx
48 real*8 :: o
49
50 do i = 1, nobs
51 o = 0.d0
52 do j = 1, dimx
53 a(i,j) = rndx()
54 b(i,j) = rndx()
55 o = o + a(i,j) * xtarg + b(i,j) * xtarg**2
56 enddo
57 obs(i) = o + noise * rndx()
58 enddo
59
60end subroutine defdata
61!
62! Main program.
63!
64!> Main program. A simple setup and call of CONOPT
65!!
66Program leastsquare
67
69 Use conopt
70 Use lsq_700
71 implicit None
72!
73! Declare the user callback routines as Integer, External:
74!
75 Integer, External :: lsq_readmatrix ! Mandatory Matrix definition routine defined below
76 Integer, External :: lsq_fdeval ! Function and Derivative evaluation routine
77 ! needed a nonlinear model.
78 Integer, External :: std_status ! Standard callback for displaying solution status
79 Integer, External :: std_solution ! Standard callback for displaying solution values
80 Integer, External :: std_message ! Standard callback for managing messages
81 Integer, External :: std_errmsg ! Standard callback for managing error messages
82#ifdef dec_directives_win32
83!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_ReadMatrix
84!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_FDEval
85!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Status
86!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Solution
87!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Message
88!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_ErrMsg
89#endif
90!
91! Control vector
92!
93 INTEGER, Dimension(:), Pointer :: cntvect
94 INTEGER :: coi_error
95
96 call startup
97!
98! Define data
99!
100 Call defdata
101!
102! Create and initialize a Control Vector
103!
104 coi_error = coi_create( cntvect )
105!
106! Tell CONOPT about the size of the model by populating the Control Vector:
107!
108 coi_error = max( coi_error, coidef_numvar( cntvect, nobs + dimx ) )
109 coi_error = max( coi_error, coidef_numcon( cntvect, nobs + 1 ) )
110 coi_error = max( coi_error, coidef_numnz( cntvect, nobs * (dimx+2) ) )
111 coi_error = max( coi_error, coidef_numnlnz( cntvect, nobs * (dimx+1) ) )
112 coi_error = max( coi_error, coidef_optdir( cntvect, -1 ) ) ! Minimize
113 coi_error = max( coi_error, coidef_objcon( cntvect, nobs + 1 ) ) ! Objective is last constraint
114 coi_error = max( coi_error, coidef_optfile( cntvect, 'leastsq.opt' ) )
115!
116! Tell CONOPT about the callback routines:
117!
118 coi_error = max( coi_error, coidef_readmatrix( cntvect, lsq_readmatrix ) )
119 coi_error = max( coi_error, coidef_fdeval( cntvect, lsq_fdeval ) )
120 coi_error = max( coi_error, coidef_status( cntvect, std_status ) )
121 coi_error = max( coi_error, coidef_solution( cntvect, std_solution ) )
122 coi_error = max( coi_error, coidef_message( cntvect, std_message ) )
123 coi_error = max( coi_error, coidef_errmsg( cntvect, std_errmsg ) )
124 coi_error = max( coi_error, coidef_debugfv( cntvect, 0 ) ) ! Debug Fdeval on or off
125
126#if defined(CONOPT_LICENSE_INT_1) && defined(CONOPT_LICENSE_INT_2) && defined(CONOPT_LICENSE_INT_3) && defined(CONOPT_LICENSE_TEXT)
127 coi_error = max( coi_error, coidef_license( cntvect, conopt_license_int_1, conopt_license_int_2, conopt_license_int_3, conopt_license_text) )
128#endif
129
130 If ( coi_error .ne. 0 ) THEN
131 write(*,*)
132 write(*,*) '**** Fatal Error while loading CONOPT Callback routines.'
133 write(*,*)
134 call flog( "Skipping Solve due to setup errors", 1 )
135 ENDIF
136!
137! Save the solution so we can check the duals:
138!
139 do_allocate = .true.
140!
141! Start CONOPT:
142!
143 coi_error = coi_solve( cntvect )
144
145 write(*,*)
146 write(*,*) 'End of Least Square example 1. Return code=',coi_error
147
148 If ( coi_error /= 0 ) then
149 call flog( "Errors encountered during solution", 1 )
150 elseif ( stacalls == 0 .or. solcalls == 0 ) then
151 call flog( "Status or Solution routine was not called", 1 )
152 elseif ( .not. ( sstat == 1 .and. mstat == 2 ) ) then
153 call flog( "Solver or Model status was not as expected (1,2)", 1 )
154 elseif ( abs( obj - 19.44434311d0 ) > 1.d-7 ) then
155 call flog( "Incorrect objective returned", 1 )
156 Else
157 Call checkdual( 'LeastSq', minimize )
158 endif
159
160 if ( coi_free(cntvect) /= 0 ) call flog( "Error while freeing control vector",1)
161
162 call flog( "Successful Solve", 0 )
163!
164! Free solution memory
165!
166 call finalize
167
168End Program leastsquare
169!
170! ============================================================================
171! Define information about the model:
172!
173
174!> Define information about the model
175!!
176!! @include{doc} readMatrix_params.dox
177Integer Function lsq_readmatrix( lower, curr, upper, vsta, type, rhs, esta, &
178 colsta, rowno, value, nlflag, n, m, nz, &
179 usrmem )
180#ifdef dec_directives_win32
181!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_ReadMatrix
182#endif
183 Use lsq_700
184 implicit none
185 integer, intent (in) :: n ! number of variables
186 integer, intent (in) :: m ! number of constraints
187 integer, intent (in) :: nz ! number of nonzeros
188 real*8, intent (in out), dimension(n) :: lower ! vector of lower bounds
189 real*8, intent (in out), dimension(n) :: curr ! vector of initial values
190 real*8, intent (in out), dimension(n) :: upper ! vector of upper bounds
191 integer, intent (in out), dimension(n) :: vsta ! vector of initial variable status
192 ! (not defined here)
193 integer, intent (out), dimension(m) :: type ! vector of equation types
194 integer, intent (in out), dimension(m) :: esta ! vector of initial equation status
195 ! (not defined here)
196 real*8, intent (in out), dimension(m) :: rhs ! vector of right hand sides
197 integer, intent (in out), dimension(n+1) :: colsta ! vector with start of column indices
198 integer, intent (out), dimension(nz) :: rowno ! vector of row numbers
199 integer, intent (in out), dimension(nz) :: nlflag ! vector of nonlinearity flags
200 real*8, intent (in out), dimension(nz) :: value ! vector of matrix values
201 real*8 usrmem(*) ! optional user memory
202
203 Integer :: i, j, k
204!
205! Information about Variables:
206! Default: Lower = -Inf, Curr = 0, and Upper = +inf.
207! Default: the status information in Vsta is not used.
208!
209 do i = 1, dimx
210 curr(i) = -0.8d0
211 enddo
212!
213! Information about Constraints:
214! Default: Rhs = 0
215! Default: the status information in Esta and the function
216! value in FV are not used.
217! Default: Type: There is no default.
218! 0 = Equality,
219! 1 = Greater than or equal,
220! 2 = Less than or equal,
221! 3 = Non binding.
222!
223! Constraints 1 to Nobs:
224! Rhs = Obs(i) and type Equality
225!
226 do i = 1, nobs
227 rhs(i) = obs(i)
228 type(i) = 0
229 enddo
230!
231! Constraint Nobs + 1 (Objective)
232! Rhs = 0 and type Non binding
233!
234 type(nobs+1) = 3
235!
236! Information about the Jacobian. CONOPT expects a columnwise
237! representation in Rowno, Value, Nlflag and Colsta.
238!
239! Colsta = Start of column indices (No Defaults):
240! Rowno = Row indices
241! Value = Value of derivative (by default only linear
242! derivatives are used)
243! Nlflag = 0 for linear and 1 for nonlinear derivative
244! (not needed for completely linear models)
245!
246!
247! Indices
248! x(i) res(j)
249! j: NL L=1
250! obj: NL
251!
252 k = 1
253 do j = 1, dimx
254 colsta(j) = k
255 do i = 1, nobs
256 rowno(k) = i
257 nlflag(k) = 1
258 k = k + 1
259 enddo
260 enddo
261 do i = 1, nobs
262 colsta(dimx+i) = k
263 rowno(k) = i
264 nlflag(k) = 0
265 value(k) = 1.d0
266 k = k + 1
267 rowno(k) = nobs+1
268 nlflag(k) = 1
269 k = k + 1
270 enddo
271 colsta(dimx+nobs+1) = k
273 lsq_readmatrix = 0 ! Return value means OK
274
275end Function lsq_readmatrix
276!
277!==========================================================================
278! Compute nonlinear terms and non-constant Jacobian elements
279!
280
281!> Compute nonlinear terms and non-constant Jacobian elements
282!!
283!! @include{doc} fdeval_params.dox
284Integer Recursive function lsq_fdeval( x, g, jac, rowno, jcnm, mode, ignerr, errcnt, &
285 n, nz, thread, usrmem )
286#ifdef dec_directives_win32
287!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_FDEval
288#endif
289 use lsq_700
290 implicit none
291 integer, intent (in) :: n ! number of variables
292 integer, intent (in) :: rowno ! number of the row to be evaluated
293 integer, intent (in) :: nz ! number of nonzeros in this row
294 real*8, intent (in), dimension(n) :: x ! vector of current solution values
295 real*8, intent (in out) :: g ! constraint value
296 real*8, intent (in out), dimension(n) :: jac ! vector of derivatives for current constraint
297 integer, intent (in), dimension(nz) :: jcnm ! list of variables that appear nonlinearly
298 ! in this row. Ffor information only.
299 integer, intent (in) :: mode ! evaluation mode: 1 = function value
300 ! 2 = derivatives, 3 = both
301 integer, intent (in) :: ignerr ! if 1 then errors can be ignored as long
302 ! as errcnt is incremented
303 integer, intent (in out) :: errcnt ! error counter to be incremented in case
304 ! of function evaluation errors.
305 integer, intent (in) :: thread
306 real*8 usrmem(*) ! optional user memory
307
308 integer :: i,j
309 real*8 :: s
310!
311! Row Nobs+1: the objective function is nonlinear
312!
313 if ( rowno .eq. nobs+1 ) then
314!
315! Mode = 1 or 3. Function value: G = sum(i, res(i)**2)
316!
317 if ( mode .eq. 1 .or. mode .eq. 3 ) then
318 s = 0.d0
319 do i = 1, nobs
320 s = s + x(dimx+i)**2
321 enddo
322 g = s
323 endif
324!
325! Mode = 2 or 3: Derivative values:
326!
327 if ( mode .eq. 2 .or. mode .eq. 3 ) then
328 do i = 1, nobs
329 jac(dimx+i) = 2*x(dimx+i)
330 enddo
331 endif
332!
333! Row 1 to Nobs: The observation definitions:
334!
335 else
336!
337! Mode = 1 or 3: Function value - x-part only
338!
339 if ( mode .eq. 1 .or. mode .eq. 3 ) then
340 s = 0.d0
341 do j = 1, dimx
342 s = s + a(rowno,j)*x(j) + b(rowno,j)*x(j)**2
343 enddo
344 g = s
345 endif
346!
347! Mode = 2 or 3: Derivatives
348!
349 if ( mode .eq. 2 .or. mode .eq. 3 ) then
350 do j = 1, dimx
351 jac(j) = a(rowno,j) + 2.d0*b(rowno,j)*x(j)
352 enddo
353 endif
354 endif
355 lsq_fdeval = 0
356
357end Function lsq_fdeval
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
subroutine checkdual(case, minmax)
Definition comdecl.f90:432
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_optfile(cntvect, optfile)
define callback routine for defining an options file.
Definition conopt.f90:928
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
#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, parameter minimize
Definition comdecl.f90:31
integer stacalls
Definition comdecl.f90:14
subroutine flog(msg, code)
Definition comdecl.f90:62
logical do_allocate
Definition comdecl.f90:27
integer mstat
Definition comdecl.f90:17
subroutine startup
Definition comdecl.f90:41