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