CONOPT
Loading...
Searching...
No Matches
leastsq5.f90
Go to the documentation of this file.
1!> @file leastsq5.f90
2!! @ingroup FORT1THREAD_EXAMPLES
3!!
4!!
5!! This model is similar to leastsq, but this time we define 2nd order
6!! information in the form of a 2DDir routine.
7!!
8!! We solve the following nonlinear least squares model:
9!!
10!! \f[
11!! \min \sum_i res_{i}^2 \\
12!!
13!! \sum_j ( a_{ij}x_j + b_{ij}x_j^2 ) + res_i = obs_i
14!! \f]
15!!
16!! where \f$a\f$, \f$b\f$, and \f$obs\f$ are known data, and \f$res\f$ and \f$x\f$ are the
17!! variables of the model.
18!!
19!!
20!! For more information about the individual callbacks, please have a look at the source code.
21
22#if defined(_WIN32) && !defined(_WIN64)
23#define dec_directives_win32
24#endif
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!!
68Program leastsquare
69
71 Use conopt
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_2ddir ! Directional 2nd derivative routine
81 Integer, External :: std_status ! Standard callback for displaying solution status
82 Integer, External :: std_solution ! Standard callback for displaying solution values
83 Integer, External :: std_message ! Standard callback for managing messages
84 Integer, External :: std_errmsg ! Standard callback for managing error messages
85#ifdef dec_directives_win32
86!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_ReadMatrix
87!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_FDEval
88!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_2DDir
89!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Status
90!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Solution
91!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Message
92!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_ErrMsg
93#endif
94!
95! Control vector
96!
97 INTEGER, Dimension(:), Pointer :: cntvect
98 INTEGER :: coi_error
99
100 call startup
101!
102! Define data
103!
104 Call defdata
105!
106! Create and initialize a Control Vector
107!
108 coi_error = coi_create( cntvect )
109!
110! Tell CONOPT about the size of the model by populating the Control Vector:
111!
112 coi_error = max( coi_error, coidef_numvar( cntvect, nobs + dimx ) )
113 coi_error = max( coi_error, coidef_numcon( cntvect, nobs + 1 ) )
114 coi_error = max( coi_error, coidef_numnz( cntvect, nobs * (dimx+2) ) )
115 coi_error = max( coi_error, coidef_numnlnz( cntvect, nobs * (dimx+1) ) )
116 coi_error = max( coi_error, coidef_optdir( cntvect, -1 ) ) ! Minimize
117 coi_error = max( coi_error, coidef_objcon( cntvect, nobs + 1 ) ) ! Objective is last constraint
118 coi_error = max( coi_error, coidef_optfile( cntvect, 'leastsq5.opt' ) )
119!
120! Tell CONOPT about the callback routines:
121!
122 coi_error = max( coi_error, coidef_readmatrix( cntvect, lsq_readmatrix ) )
123 coi_error = max( coi_error, coidef_fdeval( cntvect, lsq_fdeval ) )
124 coi_error = max( coi_error, coidef_2ddir( cntvect, lsq_2ddir ) )
125 coi_error = max( coi_error, coidef_status( cntvect, std_status ) )
126 coi_error = max( coi_error, coidef_solution( cntvect, std_solution ) )
127 coi_error = max( coi_error, coidef_message( cntvect, std_message ) )
128 coi_error = max( coi_error, coidef_errmsg( cntvect, std_errmsg ) )
129 coi_error = max( coi_error, coidef_debugfv( cntvect, 0 ) ) ! Debug Fdeval on or off
130
131#if defined(CONOPT_LICENSE_INT_1) && defined(CONOPT_LICENSE_INT_2) && defined(CONOPT_LICENSE_INT_3) && defined(CONOPT_LICENSE_TEXT)
132 coi_error = max( coi_error, coidef_license( cntvect, conopt_license_int_1, conopt_license_int_2, conopt_license_int_3, conopt_license_text) )
133#endif
134
135 If ( coi_error .ne. 0 ) THEN
136 write(*,*)
137 write(*,*) '**** Fatal Error while loading CONOPT Callback routines.'
138 write(*,*)
139 call flog( "Skipping Solve due to setup errors", 1 )
140 ENDIF
141!
142! Start CONOPT:
143!
144 coi_error = coi_solve( cntvect )
145
146 write(*,*)
147 write(*,*) 'End of Least Square example 5. Return code=',coi_error
148
149 If ( coi_error /= 0 ) then
150 call flog( "Errors encountered during solution", 1 )
151 elseif ( stacalls == 0 .or. solcalls == 0 ) then
152 call flog( "Status or Solution routine was not called", 1 )
153 elseif ( sstat /= 1 .or. mstat /= 2 ) then
154 call flog( "Solver and Model Status was not as expected (1,2)", 1 )
155 elseif ( abs( obj - 19.44434311d0 ) > 1.d-7 ) then
156 call flog( "Incorrect objective returned", 1 )
157 endif
158
159 if ( coi_free(cntvect) /= 0 ) call flog( "Error while freeing control vector",1)
160
161 call flog( "Successful Solve", 0 )
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#ifdef dec_directives_win32
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. CONOPT expects a columnwise
232! representation in Rowno, Value, Nlflag and Colsta.
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! 3:
247!
248 k = 1
249 do j = 1, dimx
250 colsta(j) = k
251 do i = 1, nobs
252 rowno(k) = i
253 nlflag(k) = 1
254 k = k + 1
255 enddo
256 enddo
257 do i = 1, nobs
258 colsta(dimx+i) = k
259 rowno(k) = i
260 nlflag(k) = 0
261 value(k) = 1.d0
262 k = k + 1
263 rowno(k) = nobs+1
264 nlflag(k) = 1
265 k = k + 1
266 enddo
267 colsta(dimx+nobs+1) = k
268
269 lsq_readmatrix = 0 ! Return value means OK
270
271end Function lsq_readmatrix
272!
273!==========================================================================
274! Compute nonlinear terms and non-constant Jacobian elements
275!
276
277!> Compute nonlinear terms and non-constant Jacobian elements
278!!
279!! @include{doc} fdeval_params.dox
280Integer Recursive function lsq_fdeval( x, g, jac, rowno, jcnm, mode, ignerr, errcnt, &
281 n, nz, thread, usrmem )
282#ifdef dec_directives_win32
283!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_FDEval
284#endif
285 use lsq_700
286 implicit none
287 integer, intent (in) :: n ! number of variables
288 integer, intent (in) :: rowno ! number of the row to be evaluated
289 integer, intent (in) :: nz ! number of nonzeros in this row
290 real*8, intent (in), dimension(n) :: x ! vector of current solution values
291 real*8, intent (in out) :: g ! constraint value
292 real*8, intent (in out), dimension(n) :: jac ! vector of derivatives for current constraint
293 integer, intent (in), dimension(nz) :: jcnm ! list of variables that appear nonlinearly
294 ! in this row. Ffor information only.
295 integer, intent (in) :: mode ! evaluation mode: 1 = function value
296 ! 2 = derivatives, 3 = both
297 integer, intent (in) :: ignerr ! if 1 then errors can be ignored as long
298 ! as errcnt is incremented
299 integer, intent (in out) :: errcnt ! error counter to be incremented in case
300 ! of function evaluation errors.
301 integer, intent (in) :: thread
302 real*8 usrmem(*) ! optional user memory
303
304 integer :: i,j
305 real*8 :: s
306!
307! Row 1: the objective function is nonlinear
308!
309 if ( rowno .eq. nobs+1 ) then
310!
311! Mode = 1 or 3. Function value: G = sum(i, res(i)**2)
312!
313 if ( mode .eq. 1 .or. mode .eq. 3 ) then
314 s = 0.d0
315 do i = 1, nobs
316 s = s + x(dimx+i)**2
317 enddo
318 g = s
319 endif
320!
321! Mode = 2 or 3: Derivative values:
322!
323 if ( mode .eq. 2 .or. mode .eq. 3 ) then
324 do i = 1, nobs
325 jac(dimx+i) = 2*x(dimx+i)
326 enddo
327 endif
328!
329! Row 1 to Nobs: The observation definitions:
330!
331 else
332!
333! Mode = 1 or 3: Function value - x-part only
334!
335 if ( mode .eq. 1 .or. mode .eq. 3 ) then
336 s = 0.d0
337 do j = 1, dimx
338 s = s + a(rowno,j)*x(j) + b(rowno,j)*x(j)**2
339 enddo
340 g = s
341 endif
342!
343! Mode = 2 or 3: Derivatives
344!
345 if ( mode .eq. 2 .or. mode .eq. 3 ) then
346 do j = 1, dimx
347 jac(j) = a(rowno,j) + 2.d0*b(rowno,j)*x(j)
348 enddo
349 endif
350 endif
351 lsq_fdeval = 0
352
353end Function lsq_fdeval
354
355
356!> Computes the second derivative of a constraint in a direction
357!!
358!! @include{doc} 2DDir_params.dox
359Integer Function lsq_2ddir( X, DX, D2G, ROWNO, JCNM, NODRV, N, NJ, THREAD, USRMEM )
360#ifdef dec_directives_win32
361!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_2DDir
362#endif
363 Use lsq_700
364 implicit none
365 INTEGER, Intent(IN) :: rowno, n, nj, thread
366 INTEGER, Intent(IN OUT) :: nodrv
367 INTEGER, Dimension(NJ), Intent(IN) :: jcnm
368 real*8, Dimension(N), Intent(IN) :: x
369 real*8, Dimension(N), Intent(IN) :: dx
370 real*8, Dimension(N), Intent(OUT) :: d2g
371 real*8, Intent(IN OUT) :: usrmem(*)
372 integer :: i, j
373
374 if ( rowno .eq. nobs+1 ) then
375!
376! Objective row: sum(i, res(i)**2 )
377!
378 do j = 1, dimx
379 d2g(j) = 0.d0
380 enddo
381 do i = 1, nobs
382 d2g(dimx+i) = 2.d0*dx(dimx+i)
383 enddo
384 else
385!
386! Constraint row with b(i,j)*x(j)**2 as only nonlinear term
387!
388 do j = 1, dimx
389 d2g(j) = 2.d0*b(rowno,j)*dx(j)
390 enddo
391 do i = 1, nobs
392 d2g(dimx+i) = 0.d0
393 enddo
394 endif
395
396 lsq_2ddir = 0
397
398End Function lsq_2ddir
integer function std_solution(xval, xmar, xbas, xsta, yval, ymar, ybas, ysta, n, m, usrmem)
Definition comdecl.f90:132
integer function std_status(modsta, solsta, iter, objval, usrmem)
Definition comdecl.f90:88
integer function std_message(smsg, dmsg, nmsg, llen, usrmem, msgv)
Definition comdecl.f90:205
integer function std_errmsg(rowno, colno, posno, msglen, usrmem, msg)
Definition comdecl.f90:248
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_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:166
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:287
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
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