CONOPT
Loading...
Searching...
No Matches
minimax.f90
Go to the documentation of this file.
1!> @file minimax.f90
2!! @ingroup FORT1THREAD_EXAMPLES
3!!
4!!
5!! We solve the following nonlinear minimax model:
6!!
7!! @verbatim
8!! min max(i, abs(res(i)) )
9!!
10!! sum(j, a(i,j)*x(j) + b(i,j)*x(j)**2 ) + res(i) = obs(i)
11!! @endverbatim
12!!
13!! where a, b, and obs are known data, and res and x are the
14!! variables of the model.
15!!
16!! To get a model with continuous derivatives we introduce an
17!! objective variable, z, and the constraints
18!! @verbatim
19!! +res(i) <= z or +res(i) - z <= 0
20!! -res(i) <= z or -res(i) - z <= 0
21!! @endverbatim
22!!
23!!
24!! For more information about the individual callbacks, please have a look at the source code.
25
26
27REAL FUNCTION rndx( )
28!
29! Defines a pseudo random number between 0 and 1
30!
31 IMPLICIT NONE
32
33 Integer, save :: seed = 12359
34
35 seed = mod(seed*1027+25,1048576)
36 rndx = float(seed)/float(1048576)
37
38END FUNCTION rndx
39
40subroutine defdata
41!
42! Define values for A, B, and Obs
43!
44 Use lsq_70
45 IMPLICIT NONE
46
47 Integer :: i, j
48 real*8, Parameter :: xtarg = -1.0
49 real*8, Parameter :: noise = 1.0
50 Real, External :: Rndx
51 real*8 :: o
52
53 do i = 1, nobs
54 o = 0.d0
55 do j = 1, dimx
56 a(i,j) = rndx()
57 b(i,j) = rndx()
58 o = o + a(i,j) * xtarg + b(i,j) * xtarg**2
59 enddo
60 obs(i) = o + noise * rndx()
61 enddo
62
63end subroutine defdata
64!
65! Main program.
66!
67!> Main program. A simple setup and call of CONOPT
68!!
69Program minimax
70
71 Use proginfo
72 Use coidef
73 Use lsq_70
74 implicit None
75!
76! Declare the user callback routines as Integer, External:
77!
78 Integer, External :: mm_readmatrix ! Mandatory Matrix definition routine defined below
79 Integer, External :: mm_fdeval ! Function and Derivative evaluation routine
80 ! needed a nonlinear model.
81 Integer, External :: std_status ! Standard callback for displaying solution status
82 Integer, External :: std_solution ! Standard callback for displaying solution values
83 Integer, External :: std_message ! Standard callback for managing messages
84 Integer, External :: std_errmsg ! Standard callback for managing error messages
85#if defined(itl)
86!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: MM_ReadMatrix
87!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: MM_FDEval
88!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Status
89!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Solution
90!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Message
91!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_ErrMsg
92#endif
93!
94! Control vector
95!
96 INTEGER :: numcallback
97 INTEGER, Dimension(:), Pointer :: cntvect
98 INTEGER :: coi_error
99
100 call startup
101!
102! Define data
103!
104 Call defdata
105!
106! Create and initialize a Control Vector
107!
108 numcallback = coidef_size()
109 Allocate( cntvect(numcallback) )
110 coi_error = coidef_inifort( cntvect )
111!
112! Tell CONOPT about the size of the model by populating the Control Vector:
113!
114 coi_error = max( coi_error, coidef_numvar( cntvect, nobs + dimx + 1 ) )
115 coi_error = max( coi_error, coidef_numcon( cntvect, 3*nobs ) )
116 coi_error = max( coi_error, coidef_numnz( cntvect, nobs * (dimx+5) ) )
117 coi_error = max( coi_error, coidef_numnlnz( cntvect, nobs * dimx ) )
118 coi_error = max( coi_error, coidef_optdir( cntvect, -1 ) ) ! Minimize
119 coi_error = max( coi_error, coidef_objvar( cntvect, nobs + dimx + 1 ) ) ! Objective is last variable
120 coi_error = max( coi_error, coidef_optfile( cntvect, 'leastsq.opt' ) )
121!
122! Tell CONOPT about the callback routines:
123!
124 coi_error = max( coi_error, coidef_readmatrix( cntvect, mm_readmatrix ) )
125 coi_error = max( coi_error, coidef_fdeval( cntvect, mm_fdeval ) )
126 coi_error = max( coi_error, coidef_status( cntvect, std_status ) )
127 coi_error = max( coi_error, coidef_solution( cntvect, std_solution ) )
128 coi_error = max( coi_error, coidef_message( cntvect, std_message ) )
129 coi_error = max( coi_error, coidef_errmsg( cntvect, std_errmsg ) )
130 coi_error = max( coi_error, coidef_debugfv( cntvect, 0 ) ) ! Debug Fdeval on or off
131
132#if defined(LICENSE_INT_1) && defined(LICENSE_INT_2) && defined(LICENSE_INT_3) && defined(LICENSE_TEXT)
133 coi_error = max( coi_error, coidef_license( cntvect, license_int_1, license_int_2, license_int_3, license_text) )
134#endif
135
136 If ( coi_error .ne. 0 ) THEN
137 write(*,*)
138 write(*,*) '**** Fatal Error while loading CONOPT Callback routines.'
139 write(*,*)
140 call flog( "Skipping Solve due to setup errors", 1 )
141 ENDIF
142!
143! Save the solution so we can check the duals:
144!
145 do_allocate = .true.
146!
147! Start CONOPT:
148!
149 coi_error = coi_solve( cntvect )
150
151 write(*,*)
152 write(*,*) 'End of Minimax example 1. Return code=',coi_error
153
154 If ( coi_error /= 0 ) then
155 call flog( "Errors encountered during solution", 1 )
156 elseif ( stacalls == 0 .or. solcalls == 0 ) then
157 call flog( "Status or Solution routine was not called", 1 )
158 elseif ( .not. ( sstat == 1 .and. mstat == 2 ) ) then
159 call flog( "Solver or Model status was not as expected (1,2)", 1 )
160! elseif ( abs( OBJ - 19.44434311d0 ) > 1.d-7 ) then
161! call flog( "Incorrect objective returned", 1 )
162 Else
163 Call checkdual( 'Minimax', minimize )
164 endif
165
166 if ( coi_free(cntvect) /= 0 ) call flog( "Error while freeing control vector",1)
167
168 call flog( "Successful Solve", 0 )
169
170End Program minimax
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 mm_readmatrix( lower, curr, upper, vsta, type, rhs, esta, &
180 colsta, rowno, value, nlflag, n, m, nz, &
181 usrmem )
182#if defined(itl)
183!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: MM_ReadMatrix
184#endif
185 Use lsq_70
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 to 3*Nobs
234! Rhs = 0 (default) and type Less than or equal
235!
236 do i = nobs+1, 3*nobs
237 type(i) = 2
238 enddo
239!
240! Information about the Jacobian. We use the standard method with
241! Rowno, Value, Nlflag and Colsta and we do not use Colno.
242!
243! Colsta = Start of column indices (No Defaults):
244! Rowno = Row indices
245! Value = Value of derivative (by default only linear
246! derivatives are used)
247! Nlflag = 0 for linear and 1 for nonlinear derivative
248! (not needed for completely linear models)
249!
250!
251! Indices
252! x(i) res(j) z
253! j: NL L=1
254! j: L=1 L=-1
255! j: L=-1 L=-1
256!
257 k = 1
258 do j = 1, dimx
259 colsta(j) = k
260 do i = 1, nobs
261 rowno(k) = i
262 nlflag(k) = 1
263 k = k + 1
264 enddo
265 enddo
266 do i = 1, nobs
267 colsta(dimx+i) = k
268 rowno(k) = i
269 nlflag(k) = 0
270 value(k) = 1.d0
271 k = k + 1
272 rowno(k) = i+nobs
273 nlflag(k) = 0
274 value(k) = 1.d0
275 k = k + 1
276 rowno(k) = i+2*nobs
277 nlflag(k) = 0
278 value(k) = -1.d0
279 k = k + 1
280 enddo
281 colsta(dimx+nobs+1) = k
282 do i = 1, 2*nobs
283 rowno(k) = nobs+i
284 nlflag(k) = 0
285 value(k) = -1.d0
286 k = k + 1
287 enddo
288 colsta(dimx+nobs+2) = k
289
290 mm_readmatrix = 0 ! Return value means OK
291
292end Function mm_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 mm_fdeval( x, g, jac, rowno, jcnm, mode, ignerr, errcnt, &
302 n, nz, thread, usrmem )
303#if defined(itl)
304!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: MM_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 nonlinear
329!
330 mm_fdeval = 0
331 if ( rowno .le. nobs ) then
332!
333! Mode = 1 or 3: Function value - x-part only
334!
335 if ( mode .eq. 1 .or. mode .eq. 3 ) then
336 s = 0.d0
337 do j = 1, dimx
338 s = s + a(rowno,j)*x(j) + b(rowno,j)*x(j)**2
339 enddo
340 g = s
341 endif
342!
343! Mode = 2 or 3: Derivatives
344!
345 if ( mode .eq. 2 .or. mode .eq. 3 ) then
346 do j = 1, dimx
347 jac(j) = a(rowno,j) + 2.d0*b(rowno,j)*x(j)
348 enddo
349 endif
350 Else
351 mm_fdeval = 1 ! Illegal row number
352 endif
353
354end 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
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
program minimax
Main program. A simple setup and call of CONOPT.
Definition minimax.f90:69
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