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