CONOPT
Loading...
Searching...
No Matches
leastl1.f90
Go to the documentation of this file.
1!> @file leastl1.f90
2!! @ingroup FORT1THREAD_EXAMPLES
3!!
4!!
5!! We solve the following nonlinear linear L1 estimation model:
6!!
7!! @verbatim
8!! min sum(i, resp(i)+resm(i) )
9!!
10!! sum(j, a(i,j)*x(j) + b(i,j)*x(j)**2 ) + resp(i) - resm(i) = obs(i)
11!! @endverbatim
12!!
13!! where a, b, and obs are known data, x are the parameter to be estimated,
14!! and resp and resm are positive variables measuring positive and
15!! negative deviations to be minimized.
16!!
17!!
18!! For more information about the individual callbacks, please have a look at the source code.
19
20
21REAL FUNCTION rndx( )
22!
23! Defines a pseudo random number between 0 and 1
24!
25 IMPLICIT NONE
26
27 Integer, save :: seed = 12359
28
29 seed = mod(seed*1027+25,1048576)
30 rndx = float(seed)/float(1048576)
31
32END FUNCTION rndx
33
34subroutine defdata
35!
36! Define values for A, B, and Obs
37!
38 Use lsq_70
39 IMPLICIT NONE
40
41 Integer :: i, j
42 real*8, Parameter :: xtarg = -1.0
43 real*8, Parameter :: noise = 1.0
44 Real, External :: Rndx
45 real*8 :: o
46
47 do i = 1, nobs
48 o = 0.d0
49 do j = 1, dimx
50 a(i,j) = rndx()
51 b(i,j) = rndx()
52 o = o + a(i,j) * xtarg + b(i,j) * xtarg**2
53 enddo
54 obs(i) = o + noise * rndx()
55 enddo
56
57end subroutine defdata
58!
59! Main program.
60!
61!> Main program. A simple setup and call of CONOPT
62!!
63Program leastl1
64
65 Use proginfo
66 Use coidef
67 Use lsq_70
68 implicit None
69!
70! Declare the user callback routines as Integer, External:
71!
72 Integer, External :: lsq_readmatrix ! Mandatory Matrix definition routine defined below
73 Integer, External :: lsq_fdeval ! Function and Derivative evaluation routine
74 ! needed a nonlinear model.
75 Integer, External :: std_status ! Standard callback for displaying solution status
76 Integer, External :: std_solution ! Standard callback for displaying solution values
77 Integer, External :: std_message ! Standard callback for managing messages
78 Integer, External :: std_errmsg ! Standard callback for managing error messages
79#if defined(itl)
80!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_ReadMatrix
81!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_FDEval
82!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Status
83!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Solution
84!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Message
85!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_ErrMsg
86#endif
87!
88! Control vector
89!
90 INTEGER :: numcallback
91 INTEGER, Dimension(:), Pointer :: cntvect
92 INTEGER :: coi_error
93
94 call startup
95!
96! Define data
97!
98 Call defdata
99!
100! Create and initialize a Control Vector
101!
102 numcallback = coidef_size()
103 Allocate( cntvect(numcallback) )
104 coi_error = coidef_inifort( cntvect )
105!
106! Tell CONOPT about the size of the model by populating the Control Vector:
107!
108 coi_error = max( coi_error, coidef_numvar( cntvect, 2*nobs + dimx ) )
109 coi_error = max( coi_error, coidef_numcon( cntvect, nobs + 1 ) )
110 coi_error = max( coi_error, coidef_numnz( cntvect, nobs * (dimx+4) ) )
111 coi_error = max( coi_error, coidef_numnlnz( cntvect, nobs * dimx ) )
112 coi_error = max( coi_error, coidef_optdir( cntvect, -1 ) ) ! Minimize
113 coi_error = max( coi_error, coidef_objcon( cntvect, nobs + 1 ) ) ! Objective is last constraint
114 coi_error = max( coi_error, coidef_optfile( cntvect, 'leastl1.opt' ) )
115!
116! Tell CONOPT about the callback routines:
117!
118 coi_error = max( coi_error, coidef_readmatrix( cntvect, lsq_readmatrix ) )
119 coi_error = max( coi_error, coidef_fdeval( cntvect, lsq_fdeval ) )
120 coi_error = max( coi_error, coidef_status( cntvect, std_status ) )
121 coi_error = max( coi_error, coidef_solution( cntvect, std_solution ) )
122 coi_error = max( coi_error, coidef_message( cntvect, std_message ) )
123 coi_error = max( coi_error, coidef_errmsg( cntvect, std_errmsg ) )
124 coi_error = max( coi_error, coidef_debugfv( cntvect, 0 ) ) ! Debug Fdeval on or off
125
126#if defined(LICENSE_INT_1) && defined(LICENSE_INT_2) && defined(LICENSE_INT_3) && defined(LICENSE_TEXT)
127 coi_error = max( coi_error, coidef_license( cntvect, license_int_1, license_int_2, license_int_3, license_text) )
128#endif
129
130 If ( coi_error .ne. 0 ) THEN
131 write(*,*)
132 write(*,*) '**** Fatal Error while loading CONOPT Callback routines.'
133 write(*,*)
134 call flog( "Skipping Solve due to setup errors", 1 )
135 ENDIF
136!
137! Save the solution so we can check the duals:
138!
139 do_allocate = .true.
140!
141! Start CONOPT:
142!
143 coi_error = coi_solve( cntvect )
144
145 write(*,*)
146 write(*,*) 'End of Least Square example 1. Return code=',coi_error
147
148 If ( coi_error /= 0 ) then
149 call flog( "Errors encountered during solution", 1 )
150 elseif ( stacalls == 0 .or. solcalls == 0 ) then
151 call flog( "Status or Solution routine was not called", 1 )
152 elseif ( .not. ( sstat == 1 .and. mstat == 2 ) ) then
153 call flog( "Solver or Model status was not as expected (1,2)", 1 )
154! elseif ( abs( OBJ - 19.44434311d0 ) > 1.d-7 ) then
155! call flog( "Incorrect objective returned", 1 )
156 Else
157 Call checkdual( 'LeastL1', minimize )
158 endif
159
160 if ( coi_free(cntvect) /= 0 ) call flog( "Error while freeing control vector",1)
161
162 call flog( "Successful Solve", 0 )
163
164End Program leastl1
165!
166! ============================================================================
167! Define information about the model:
168!
169
170!> Define information about the model
171!!
172!! @include{doc} readMatrix_params.dox
173Integer Function lsq_readmatrix( lower, curr, upper, vsta, type, rhs, esta, &
174 colsta, rowno, value, nlflag, n, m, nz, &
175 usrmem )
176#if defined(itl)
177!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_ReadMatrix
178#endif
179 Use lsq_70
180 implicit none
181 integer, intent (in) :: n ! number of variables
182 integer, intent (in) :: m ! number of constraints
183 integer, intent (in) :: nz ! number of nonzeros
184 real*8, intent (in out), dimension(n) :: lower ! vector of lower bounds
185 real*8, intent (in out), dimension(n) :: curr ! vector of initial values
186 real*8, intent (in out), dimension(n) :: upper ! vector of upper bounds
187 integer, intent (in out), dimension(n) :: vsta ! vector of initial variable status
188 ! (not defined here)
189 integer, intent (out), dimension(m) :: type ! vector of equation types
190 integer, intent (in out), dimension(m) :: esta ! vector of initial equation status
191 ! (not defined here)
192 real*8, intent (in out), dimension(m) :: rhs ! vector of right hand sides
193 integer, intent (in out), dimension(n+1) :: colsta ! vector with start of column indices
194 integer, intent (out), dimension(nz) :: rowno ! vector of row numbers
195 integer, intent (in out), dimension(nz) :: nlflag ! vector of nonlinearity flags
196 real*8, intent (in out), dimension(nz) :: value ! vector of matrix values
197 real*8 usrmem(*) ! optional user memory
198
199 Integer :: i, j, k
200!
201! Information about Variables:
202! Default: Lower = -Inf, Curr = 0, and Upper = +inf.
203! Default: the status information in Vsta is not used.
204! X : New initial = 0.8
205! Resp: New lower = 0.0
206! Resm: New lower = 0.0
207!
208 do i = 1, dimx
209 curr(i) = -0.8d0
210 enddo
211 do j = 1, 2*nobs
212 lower(dimx+j) = 0.0d0
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. We use the standard method with
239! Rowno, Value, Nlflag and Colsta and we do not use Colno.
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
286
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#if defined(itl)
301!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_FDEval
302#endif
303 use lsq_70
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 :: 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: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
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