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