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