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