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