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!! \f[
8!! \min \sum_i ( resp_i + resm_i ) \\
9!! \sum_j (a_{ij}x_j + b_{ij}x_{j}^2 ) + resp_{i} - resm_{i} = obs_{i}
10!! \f]
11!!
12!! where \f$a\f$, \f$b\f$, and \f$obs\f$ are known data, \f$x\f$ are the parameter to be estimated,
13!! and \f$resp\f$ and \f$resm\f$ are positive variables measuring positive and
14!! negative deviations to be minimized.
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
23REAL FUNCTION rndx( )
24!
25! Defines a pseudo random number between 0 and 1
26!
27 IMPLICIT NONE
28
29 Integer, save :: seed = 12359
30
31 seed = mod(seed*1027+25,1048576)
32 rndx = float(seed)/float(1048576)
33
34END FUNCTION rndx
35
36subroutine defdata
37!
38! Define values for A, B, and Obs
39!
40 Use lsq_70
41 IMPLICIT NONE
42
43 Integer :: i, j
44 real*8, Parameter :: xtarg = -1.0
45 real*8, Parameter :: noise = 1.0
46 Real, External :: Rndx
47 real*8 :: o
48
49 do i = 1, nobs
50 o = 0.d0
51 do j = 1, dimx
52 a(i,j) = rndx()
53 b(i,j) = rndx()
54 o = o + a(i,j) * xtarg + b(i,j) * xtarg**2
55 enddo
56 obs(i) = o + noise * rndx()
57 enddo
58
59end subroutine defdata
60!
61! Main program.
62!
63!> Main program. A simple setup and call of CONOPT
64!!
65Program leastl1
66
68 Use conopt
69 Use lsq_70
70 implicit None
71!
72! Declare the user callback routines as Integer, External:
73!
74 Integer, External :: lsq_readmatrix ! Mandatory Matrix definition routine defined below
75 Integer, External :: lsq_fdeval ! Function and Derivative evaluation routine
76 ! needed a nonlinear model.
77 Integer, External :: std_status ! Standard callback for displaying solution status
78 Integer, External :: std_solution ! Standard callback for displaying solution values
79 Integer, External :: std_message ! Standard callback for managing messages
80 Integer, External :: std_errmsg ! Standard callback for managing error messages
81#ifdef dec_directives_win32
82!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_ReadMatrix
83!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_FDEval
84!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Status
85!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Solution
86!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Message
87!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_ErrMsg
88#endif
89!
90! Control vector
91!
92 INTEGER, Dimension(:), Pointer :: cntvect
93 INTEGER :: coi_error
94
95 call startup
96!
97! Define data
98!
99 Call defdata
100!
101! Create and initialize a Control Vector
102!
103 coi_error = coi_create( cntvect )
104!
105! Tell CONOPT about the size of the model by populating the Control Vector:
106!
107 coi_error = max( coi_error, coidef_numvar( cntvect, 2*nobs + dimx ) )
108 coi_error = max( coi_error, coidef_numcon( cntvect, nobs + 1 ) )
109 coi_error = max( coi_error, coidef_numnz( cntvect, nobs * (dimx+4) ) )
110 coi_error = max( coi_error, coidef_numnlnz( cntvect, nobs * dimx ) )
111 coi_error = max( coi_error, coidef_optdir( cntvect, -1 ) ) ! Minimize
112 coi_error = max( coi_error, coidef_objcon( cntvect, nobs + 1 ) ) ! Objective is last constraint
113 coi_error = max( coi_error, coidef_optfile( cntvect, 'leastl1.opt' ) )
114!
115! Tell CONOPT about the callback routines:
116!
117 coi_error = max( coi_error, coidef_readmatrix( cntvect, lsq_readmatrix ) )
118 coi_error = max( coi_error, coidef_fdeval( cntvect, lsq_fdeval ) )
119 coi_error = max( coi_error, coidef_status( cntvect, std_status ) )
120 coi_error = max( coi_error, coidef_solution( cntvect, std_solution ) )
121 coi_error = max( coi_error, coidef_message( cntvect, std_message ) )
122 coi_error = max( coi_error, coidef_errmsg( cntvect, std_errmsg ) )
123 coi_error = max( coi_error, coidef_debugfv( cntvect, 0 ) ) ! Debug Fdeval on or off
124
125#if defined(CONOPT_LICENSE_INT_1) && defined(CONOPT_LICENSE_INT_2) && defined(CONOPT_LICENSE_INT_3) && defined(CONOPT_LICENSE_TEXT)
126 coi_error = max( coi_error, coidef_license( cntvect, conopt_license_int_1, conopt_license_int_2, conopt_license_int_3, conopt_license_text) )
127#endif
128
129 If ( coi_error .ne. 0 ) THEN
130 write(*,*)
131 write(*,*) '**** Fatal Error while loading CONOPT Callback routines.'
132 write(*,*)
133 call flog( "Skipping Solve due to setup errors", 1 )
134 ENDIF
135!
136! Save the solution so we can check the duals:
137!
138 do_allocate = .true.
139!
140! Start CONOPT:
141!
142 coi_error = coi_solve( cntvect )
143
144 write(*,*)
145 write(*,*) 'End of Least Square example 1. Return code=',coi_error
146
147 If ( coi_error /= 0 ) then
148 call flog( "Errors encountered during solution", 1 )
149 elseif ( stacalls == 0 .or. solcalls == 0 ) then
150 call flog( "Status or Solution routine was not called", 1 )
151 elseif ( .not. ( sstat == 1 .and. mstat == 2 ) ) then
152 call flog( "Solver or Model status was not as expected (1,2)", 1 )
153! elseif ( abs( OBJ - 19.44434311d0 ) > 1.d-7 ) then
154! call flog( "Incorrect objective returned", 1 )
155 Else
156 Call checkdual( 'LeastL1', minimize )
157 endif
158
159 if ( coi_free(cntvect) /= 0 ) call flog( "Error while freeing control vector",1)
160
161 call flog( "Successful Solve", 0 )
162!
163! Free solution memory
164!
165 call finalize
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#ifdef dec_directives_win32
180!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_ReadMatrix
181#endif
182 Use lsq_70
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! X : New initial = 0.8
208! Resp: New lower = 0.0
209! Resm: New lower = 0.0
210!
211 do i = 1, dimx
212 curr(i) = -0.8d0
213 enddo
214 do j = 1, 2*nobs
215 lower(dimx+j) = 0.0d0
216 Enddo
217!
218! Information about Constraints:
219! Default: Rhs = 0
220! Default: the status information in Esta and the function
221! value in FV are not used.
222! Default: Type: There is no default.
223! 0 = Equality,
224! 1 = Greater than or equal,
225! 2 = Less than or equal,
226! 3 = Non binding.
227!
228! Constraints 1 to Nobs:
229! Rhs = Obs(i) and type Equality
230!
231 do i = 1, nobs
232 rhs(i) = obs(i)
233 type(i) = 0
234 enddo
235!
236! Constraint Nobs + 1 (Objective)
237! Rhs = 0 and type Non binding
238!
239 type(nobs+1) = 3
240!
241! Information about the Jacobian. CONOPT expects a columnwise
242! representation in Rowno, Value, Nlflag and Colsta.
243!
244! Colsta = Start of column indices (No Defaults):
245! Rowno = Row indices
246! Value = Value of derivative (by default only linear
247! derivatives are used)
248! Nlflag = 0 for linear and 1 for nonlinear derivative
249! (not needed for completely linear models)
250!
251!
252! Indices
253! x(i) resp(j) resm(j)
254! j: NL L=1 L=-1
255! obj: L=1 L=1
256!
257 k = 1
258 do i = 1, dimx ! x(i) variables
259 colsta(i) = k
260 do j = 1, nobs
261 rowno(k) = j
262 nlflag(k) = 1
263 k = k + 1
264 enddo
265 enddo
266 do j = 1, nobs ! resp(j) variables
267 colsta(dimx+j) = k
268 rowno(k) = j
269 nlflag(k) = 0
270 value(k) = 1.d0
271 k = k + 1
272 rowno(k) = nobs+1
273 nlflag(k) = 0
274 value(k) = 1.d0
275 k = k + 1
276 enddo
277 do j = 1, nobs ! resm(j) variables
278 colsta(dimx+nobs+j) = k
279 rowno(k) = j
280 nlflag(k) = 0
281 value(k) = -1.d0
282 k = k + 1
283 rowno(k) = nobs+1
284 nlflag(k) = 0
285 value(k) = 1.d0
286 k = k + 1
287 enddo
288 colsta(dimx+2*nobs+1) = k
290 lsq_readmatrix = 0 ! Return value means OK
291
292end Function lsq_readmatrix
293!
294!==========================================================================
295! Compute nonlinear terms and non-constant Jacobian elements
296!
297
298!> Compute nonlinear terms and non-constant Jacobian elements
299!!
300!! @include{doc} fdeval_params.dox
301Integer Function lsq_fdeval( x, g, jac, rowno, jcnm, mode, ignerr, errcnt, &
302 n, nz, thread, usrmem )
303#ifdef dec_directives_win32
304!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_FDEval
305#endif
306 use lsq_70
307 implicit none
308 integer, intent (in) :: n ! number of variables
309 integer, intent (in) :: rowno ! number of the row to be evaluated
310 integer, intent (in) :: nz ! number of nonzeros in this row
311 real*8, intent (in), dimension(n) :: x ! vector of current solution values
312 real*8, intent (in out) :: g ! constraint value
313 real*8, intent (in out), dimension(n) :: jac ! vector of derivatives for current constraint
314 integer, intent (in), dimension(nz) :: jcnm ! list of variables that appear nonlinearly
315 ! in this row. Ffor information only.
316 integer, intent (in) :: mode ! evaluation mode: 1 = function value
317 ! 2 = derivatives, 3 = both
318 integer, intent (in) :: ignerr ! if 1 then errors can be ignored as long
319 ! as errcnt is incremented
320 integer, intent (in out) :: errcnt ! error counter to be incremented in case
321 ! of function evaluation errors.
322 integer, intent (in) :: thread
323 real*8 usrmem(*) ! optional user memory
324
325 integer :: j
326 real*8 :: s
327!
328! Row Nobs+1: the objective function is linear
329!
330 if ( rowno .eq. nobs+1 ) then
331 lsq_fdeval = 1; Return ! This is an error
332!
333! Row 1 to Nobs: The observation definitions:
334!
335 else
336!
337! Mode = 1 or 3: Function value - x-part only
338!
339 if ( mode .eq. 1 .or. mode .eq. 3 ) then
340 s = 0.d0
341 do j = 1, dimx
342 s = s + a(rowno,j)*x(j) + b(rowno,j)*x(j)**2
343 enddo
344 g = s
345 endif
346!
347! Mode = 2 or 3: Derivatives
348!
349 if ( mode .eq. 2 .or. mode .eq. 3 ) then
350 do j = 1, dimx
351 jac(j) = a(rowno,j) + 2.d0*b(rowno,j)*x(j)
352 enddo
353 endif
354 endif
355 lsq_fdeval = 0
356
357end 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
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