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