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