CONOPT
Loading...
Searching...
No Matches
minimax01.f90
Go to the documentation of this file.
1!> @file minimax01.f90
2!! @ingroup FORT1THREAD_EXAMPLES
3!!
4!!
5!! MiniMax01 is a small version of minimax used for quick turnaround.
6!!
7!! We solve the following nonlinear minimax model:
8!!
9!! @verbatim
10!! min max(i, abs(res(i)) )
11!!
12!! sum(j, a(i,j)*x(j) + b(i,j)*x(j)**2 ) + res(i) = obs(i)
13!! @endverbatim
14!!
15!! where a, b, and obs are known data, and res and x are the
16!! variables of the model.
17!!
18!! To get a model with continuous derivatives we introduce an
19!! objective variable, z, and the constraints
20!! @verbatim
21!! +res(i) <= z or +res(i) - z <= 0
22!! -res(i) <= z or -res(i) - z <= 0
23!! @endverbatim
24!!
25!!
26!! For more information about the individual callbacks, please have a look at the source code.
27
28
29REAL FUNCTION rndx( )
30!
31! Defines a pseudo random number between 0 and 1
32!
33 IMPLICIT NONE
34
35 Integer, save :: seed = 12359
36
37 seed = mod(seed*1027+25,1048576)
38 rndx = float(seed)/float(1048576)
39
40END FUNCTION rndx
41
42subroutine defdata
43!
44! Define values for A, B, and Obs
45!
46 Use lsq_7
47 IMPLICIT NONE
48
49 Integer :: i, j
50 real*8, Parameter :: xtarg = -1.0
51 real*8, Parameter :: noise = 1.0
52 Real, External :: Rndx
53 real*8 :: o
54
55 do i = 1, nobs
56 o = 0.d0
57 do j = 1, dimx
58 a(i,j) = rndx()
59 b(i,j) = rndx()
60 o = o + a(i,j) * xtarg + b(i,j) * xtarg**2
61 enddo
62 obs(i) = o + noise * rndx()
63 enddo
64
65end subroutine defdata
66!
67! Main program.
68!
69!> Main program. A simple setup and call of CONOPT
70!!
71Program minimax01
72
73 Use proginfo
74 Use coidef
75 Use lsq_7
76 implicit None
77!
78! Declare the user callback routines as Integer, External:
79!
80 Integer, External :: mm_readmatrix ! Mandatory Matrix definition routine defined below
81 Integer, External :: mm_fdeval ! Function and Derivative evaluation routine
82 ! needed a nonlinear model.
83 Integer, External :: std_status ! Standard callback for displaying solution status
84 Integer, External :: std_solution ! Standard callback for displaying solution values
85 Integer, External :: std_message ! Standard callback for managing messages
86 Integer, External :: std_errmsg ! Standard callback for managing error messages
87#if defined(itl)
88!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: MM_ReadMatrix
89!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: MM_FDEval
90!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Status
91!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Solution
92!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Message
93!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_ErrMsg
94#endif
95!
96! Control vector
97!
98 INTEGER :: numcallback
99 INTEGER, Dimension(:), Pointer :: cntvect
100 INTEGER :: coi_error
101
102 call startup
103!
104! Define data
105!
106 Call defdata
107!
108! Create and initialize a Control Vector
109!
110 numcallback = coidef_size()
111 Allocate( cntvect(numcallback) )
112 coi_error = coidef_inifort( cntvect )
113!
114! Tell CONOPT about the size of the model by populating the Control Vector:
115!
116 coi_error = max( coi_error, coidef_numvar( cntvect, nobs + dimx + 1 ) )
117 coi_error = max( coi_error, coidef_numcon( cntvect, 3*nobs ) )
118 coi_error = max( coi_error, coidef_numnz( cntvect, nobs * (dimx+5) ) )
119 coi_error = max( coi_error, coidef_numnlnz( cntvect, nobs * dimx ) )
120 coi_error = max( coi_error, coidef_optdir( cntvect, -1 ) ) ! Minimize
121 coi_error = max( coi_error, coidef_objvar( cntvect, nobs + dimx + 1 ) ) ! Objective is last variable
122 coi_error = max( coi_error, coidef_optfile( cntvect, 'MiniMax01.opt' ) )
123!
124! Tell CONOPT about the callback routines:
125!
126 coi_error = max( coi_error, coidef_readmatrix( cntvect, mm_readmatrix ) )
127 coi_error = max( coi_error, coidef_fdeval( cntvect, mm_fdeval ) )
128 coi_error = max( coi_error, coidef_status( cntvect, std_status ) )
129 coi_error = max( coi_error, coidef_solution( cntvect, std_solution ) )
130 coi_error = max( coi_error, coidef_message( cntvect, std_message ) )
131 coi_error = max( coi_error, coidef_errmsg( cntvect, std_errmsg ) )
132 coi_error = max( coi_error, coidef_debugfv( cntvect, 0 ) ) ! Debug Fdeval on or off
133
134#if defined(LICENSE_INT_1) && defined(LICENSE_INT_2) && defined(LICENSE_INT_3) && defined(LICENSE_TEXT)
135 coi_error = max( coi_error, coidef_license( cntvect, license_int_1, license_int_2, license_int_3, license_text) )
136#endif
137
138 If ( coi_error .ne. 0 ) THEN
139 write(*,*)
140 write(*,*) '**** Fatal Error while loading CONOPT Callback routines.'
141 write(*,*)
142 call flog( "Skipping Solve due to setup errors", 1 )
143 ENDIF
144!
145! Save the solution so we can check the duals:
146!
147 do_allocate = .true.
148!
149! Start CONOPT:
150!
151 coi_error = coi_solve( cntvect )
152
153 write(*,*)
154 write(*,*) 'End of MiniMax01 example 1. Return code=',coi_error
155
156 If ( coi_error /= 0 ) then
157 call flog( "Errors encountered during solution", 1 )
158 elseif ( stacalls == 0 .or. solcalls == 0 ) then
159 call flog( "Status or Solution routine was not called", 1 )
160 elseif ( .not. ( sstat == 1 .and. mstat == 2 ) ) then
161 call flog( "Solver or Model status was not as expected (1,2)", 1 )
162! elseif ( abs( OBJ - 19.44434311d0 ) > 1.d-7 ) then
163! call flog( "Incorrect objective returned", 1 )
164 Else
165 Call checkdual( 'MiniMax01', minimize )
166 endif
167
168 if ( coi_free(cntvect) /= 0 ) call flog( "Error while freeing control vector",1)
169
170 call flog( "Successful Solve", 0 )
171
172End Program minimax01
173!
174! ============================================================================
175! Define information about the model:
176!
177
178!> Define information about the model
179!!
180!! @include{doc} readMatrix_params.dox
181Integer Function mm_readmatrix( lower, curr, upper, vsta, type, rhs, esta, &
182 colsta, rowno, value, nlflag, n, m, nz, &
183 usrmem )
184#if defined(itl)
185!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: MM_ReadMatrix
186#endif
187 Use lsq_7
188 implicit none
189 integer, intent (in) :: n ! number of variables
190 integer, intent (in) :: m ! number of constraints
191 integer, intent (in) :: nz ! number of nonzeros
192 real*8, intent (in out), dimension(n) :: lower ! vector of lower bounds
193 real*8, intent (in out), dimension(n) :: curr ! vector of initial values
194 real*8, intent (in out), dimension(n) :: upper ! vector of upper bounds
195 integer, intent (in out), dimension(n) :: vsta ! vector of initial variable status
196 ! (not defined here)
197 integer, intent (out), dimension(m) :: type ! vector of equation types
198 integer, intent (in out), dimension(m) :: esta ! vector of initial equation status
199 ! (not defined here)
200 real*8, intent (in out), dimension(m) :: rhs ! vector of right hand sides
201 integer, intent (in out), dimension(n+1) :: colsta ! vector with start of column indices
202 integer, intent (out), dimension(nz) :: rowno ! vector of row numbers
203 integer, intent (in out), dimension(nz) :: nlflag ! vector of nonlinearity flags
204 real*8, intent (in out), dimension(nz) :: value ! vector of matrix values
205 real*8 usrmem(*) ! optional user memory
206
207 Integer :: i, j, k
208!
209! Information about Variables:
210! Default: Lower = -Inf, Curr = 0, and Upper = +inf.
211! Default: the status information in Vsta is not used.
212!
213 do j = 1, dimx
214 curr(j) = -0.8d0
215 enddo
216!
217! Information about Constraints:
218! Default: Rhs = 0
219! Default: the status information in Esta and the function
220! value in FV are not used.
221! Default: Type: There is no default.
222! 0 = Equality,
223! 1 = Greater than or equal,
224! 2 = Less than or equal,
225! 3 = Non binding.
226!
227! Constraints 1 to Nobs:
228! Rhs = Obs(i) and type Equality
229!
230 do i = 1, nobs
231 rhs(i) = obs(i)
232 type(i) = 0
233 enddo
234!
235! Constraint Nobs+1 to 3*Nobs
236! Rhs = 0 (default) and type Less than or equal
237!
238 do i = nobs+1, 3*nobs
239 type(i) = 2
240 enddo
241!
242! Information about the Jacobian. We use the standard method with
243! Rowno, Value, Nlflag and Colsta and we do not use Colno.
244!
245! Colsta = Start of column indices (No Defaults):
246! Rowno = Row indices
247! Value = Value of derivative (by default only linear
248! derivatives are used)
249! Nlflag = 0 for linear and 1 for nonlinear derivative
250! (not needed for completely linear models)
251!
252!
253! Indices
254! x(j) res(i) z
255! i: NL L=1
256! i: L=1 L=-1
257! i: L=-1 L=-1
258!
259 k = 1
260 do j = 1, dimx
261 colsta(j) = k
262 do i = 1, nobs
263 rowno(k) = i
264 nlflag(k) = 1
265 k = k + 1
266 enddo
267 enddo
268 do i = 1, nobs
269 colsta(dimx+i) = k
270 rowno(k) = i
271 nlflag(k) = 0
272 value(k) = 1.d0
273 k = k + 1
274 rowno(k) = i+nobs
275 nlflag(k) = 0
276 value(k) = 1.d0
277 k = k + 1
278 rowno(k) = i+2*nobs
279 nlflag(k) = 0
280 value(k) = -1.d0
281 k = k + 1
282 enddo
283 colsta(dimx+nobs+1) = k
284 do i = 1, 2*nobs
285 rowno(k) = nobs+i
286 nlflag(k) = 0
287 value(k) = -1.d0
288 k = k + 1
289 enddo
290 colsta(dimx+nobs+2) = k
291
292 mm_readmatrix = 0 ! Return value means OK
293
294end Function mm_readmatrix
295!
296!==========================================================================
297! Compute nonlinear terms and non-constant Jacobian elements
298!
299
300!> Compute nonlinear terms and non-constant Jacobian elements
301!!
302!! @include{doc} fdeval_params.dox
303Integer Function mm_fdeval( x, g, jac, rowno, jcnm, mode, ignerr, errcnt, &
304 n, nz, thread, usrmem )
305#if defined(itl)
306!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: MM_FDEval
307#endif
308 use lsq_7
309 implicit none
310 integer, intent (in) :: n ! number of variables
311 integer, intent (in) :: rowno ! number of the row to be evaluated
312 integer, intent (in) :: nz ! number of nonzeros in this row
313 real*8, intent (in), dimension(n) :: x ! vector of current solution values
314 real*8, intent (in out) :: g ! constraint value
315 real*8, intent (in out), dimension(n) :: jac ! vector of derivatives for current constraint
316 integer, intent (in), dimension(nz) :: jcnm ! list of variables that appear nonlinearly
317 ! in this row. Ffor information only.
318 integer, intent (in) :: mode ! evaluation mode: 1 = function value
319 ! 2 = derivatives, 3 = both
320 integer, intent (in) :: ignerr ! if 1 then errors can be ignored as long
321 ! as errcnt is incremented
322 integer, intent (in out) :: errcnt ! error counter to be incremented in case
323 ! of function evaluation errors.
324 integer, intent (in) :: thread
325 real*8 usrmem(*) ! optional user memory
326
327 integer :: j
328 real*8 :: s
329!
330! Row Nobs+1: the objective function is nonlinear
331!
332 mm_fdeval = 0
333 if ( rowno .le. nobs ) then
334!
335! Mode = 1 or 3: Function value - x-part only
336!
337 if ( mode .eq. 1 .or. mode .eq. 3 ) then
338 s = 0.d0
339 do j = 1, dimx
340 s = s + a(rowno,j)*x(j) + b(rowno,j)*x(j)**2
341 enddo
342 g = s
343 endif
344!
345! Mode = 2 or 3: Derivatives
346!
347 if ( mode .eq. 2 .or. mode .eq. 3 ) then
348 do j = 1, dimx
349 jac(j) = a(rowno,j) + 2.d0*b(rowno,j)*x(j)
350 enddo
351 endif
352 Else
353 mm_fdeval = 1 ! Illegal row number
354 endif
355
356end Function mm_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_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_objvar(cntvect, objvar)
defines the Objective Variable.
Definition coistart.f90:586
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
#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 minimax01
Main program. A simple setup and call of CONOPT.
Definition minimax01.f90:71
integer function mm_fdeval(x, g, jac, rowno, jcnm, mode, ignerr, errcnt, n, nz, thread, usrmem)
Compute nonlinear terms and non-constant Jacobian elements.
Definition minimax.f90:303
integer function mm_readmatrix(lower, curr, upper, vsta, type, rhs, esta, colsta, rowno, value, nlflag, n, m, nz, usrmem)
Define information about the model.
Definition minimax.f90:182
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