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