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