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!! @verbatim
11!! k=1: z(1)
12!! k=2: z(1)**2 + z(2)
13!! k=3: z(2)**2 + z(3)
14!! @endverbatim
15!!
16!! We solve the following nonlinear minimax model:
17!!
18!! @verbatim
19!! min sum(k, max(i, abs(res(i,k)) )
20!!
21!! sum(j, a(i,j,k)*x(j,k) + b(i,j,k)*x(j,k)**2 ) + res(i,k) = obs(i,k)
22!! all k: sum(j,x(j,k)) =L= 1;
23!! @endverbatim
24!!
25!! where a, b, and obs are known data, and res and x are the
26!! variables of the model.
27!!
28!! To get a model with continuous derivatives we introduce an
29!! objective variable, z, and the constraints
30!! @verbatim
31!! +res(i,k) <= z(k) or +res(i,k) - z(k) <= 0
32!! -res(i,k) <= z(k) or +res(i,k) + z(k) >= 0
33!! @endverbatim
34!! Note opposite sign from Minimax01
35!! and minimize sum(k, z(k) )
36!!
37!!
38!! For more information about the individual callbacks, please have a look at the source code.
39
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 minimax07
86
87 Use proginfo
88 Use coidef
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#if defined(itl)
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 :: numcallback
113 INTEGER, Dimension(:), Pointer :: cntvect
114 INTEGER :: coi_error
115
116 call startup
117!
118! Define data
119!
120 Call defdata
121!
122! Create and initialize a Control Vector
123!
124 numcallback = coidef_size()
125 Allocate( cntvect(numcallback) )
126 coi_error = coidef_inifort( cntvect )
127!
128! Tell CONOPT about the size of the model by populating the Control Vector:
129!
130 coi_error = max( coi_error, coidef_numvar( cntvect, nobs*dimk + dimx*dimk + dimk ) )
131 coi_error = max( coi_error, coidef_numcon( cntvect, 3*nobs*dimk+dimk+1 ) )
132 coi_error = max( coi_error, coidef_numnz( cntvect, nobs * dimx * dimk + 5* nobs*dimk + dimk + dimk*dimx + 4*nobs ) )
133 coi_error = max( coi_error, coidef_numnlnz( cntvect, nobs * dimx * dimk + 4*nobs ) )
134 coi_error = max( coi_error, coidef_optdir( cntvect, -1 ) ) ! Minimize
135 coi_error = max( coi_error, coidef_objcon( cntvect, 3*nobs*dimk + dimk+ 1 ) ) ! Objective is last constraint
136 coi_error = max( coi_error, coidef_optfile( cntvect, 'Minimax07.opt' ) )
137!
138! Tell CONOPT about the callback routines:
139!
140 coi_error = max( coi_error, coidef_readmatrix( cntvect, mm_readmatrix ) )
141 coi_error = max( coi_error, coidef_fdeval( cntvect, mm_fdeval ) )
142 coi_error = max( coi_error, coidef_status( cntvect, std_status ) )
143 coi_error = max( coi_error, coidef_solution( cntvect, std_solution ) )
144 coi_error = max( coi_error, coidef_message( cntvect, std_message ) )
145 coi_error = max( coi_error, coidef_errmsg( cntvect, std_errmsg ) )
146 coi_error = max( coi_error, coidef_debugfv( cntvect, 0 ) ) ! Debug Fdeval on or off
147
148#if defined(LICENSE_INT_1) && defined(LICENSE_INT_2) && defined(LICENSE_INT_3) && defined(LICENSE_TEXT)
149 coi_error = max( coi_error, coidef_license( cntvect, license_int_1, license_int_2, license_int_3, license_text) )
150#endif
151
152 If ( coi_error .ne. 0 ) THEN
153 write(*,*)
154 write(*,*) '**** Fatal Error while loading CONOPT Callback routines.'
155 write(*,*)
156 call flog( "Skipping Solve due to setup errors", 1 )
157 ENDIF
158!
159! Save the solution so we can check the duals:
160!
161 do_allocate = .true.
162!
163! Start CONOPT:
164!
165 coi_error = coi_solve( cntvect )
166
167 write(*,*)
168 write(*,*) 'End of Minimax07 example 1. Return code=',coi_error
169
170 If ( coi_error /= 0 ) then
171 call flog( "Errors encountered during solution", 1 )
172 elseif ( stacalls == 0 .or. solcalls == 0 ) then
173 call flog( "Status or Solution routine was not called", 1 )
174 elseif ( .not. ( sstat == 1 .and. mstat == 2 ) ) then
175 call flog( "Solver or Model status was not as expected (1,2)", 1 )
176! elseif ( abs( OBJ - 19.44434311d0 ) > 1.d-7 ) then
177! call flog( "Incorrect objective returned", 1 )
178 Else
179 Call checkdual( 'Minimax07', minimize )
180 endif
181
182 if ( coi_free(cntvect) /= 0 ) call flog( "Error while freeing control vector",1)
183
184 call flog( "Successful Solve", 0 )
185
186End Program minimax07
187!
188! ============================================================================
189! Define information about the model:
190!
191
192!> Define information about the model
193!!
194!! @include{doc} readMatrix_params.dox
195Integer Function mm_readmatrix( lower, curr, upper, vsta, type, rhs, esta, &
196 colsta, rowno, value, nlflag, n, m, nz, &
197 usrmem )
198#if defined(itl)
199!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: MM_ReadMatrix
200#endif
201 Use lsq_7k
202 implicit none
203 integer, intent (in) :: n ! number of variables
204 integer, intent (in) :: m ! number of constraints
205 integer, intent (in) :: nz ! number of nonzeros
206 real*8, intent (in out), dimension(n) :: lower ! vector of lower bounds
207 real*8, intent (in out), dimension(n) :: curr ! vector of initial values
208 real*8, intent (in out), dimension(n) :: upper ! vector of upper bounds
209 integer, intent (in out), dimension(n) :: vsta ! vector of initial variable status
210 ! (not defined here)
211 integer, intent (out), dimension(m) :: type ! vector of equation types
212 integer, intent (in out), dimension(m) :: esta ! vector of initial equation status
213 ! (not defined here)
214 real*8, intent (in out), dimension(m) :: rhs ! vector of right hand sides
215 integer, intent (in out), dimension(n+1) :: colsta ! vector with start of column indices
216 integer, intent (out), dimension(nz) :: rowno ! vector of row numbers
217 integer, intent (in out), dimension(nz) :: nlflag ! vector of nonlinearity flags
218 real*8, intent (in out), dimension(nz) :: value ! vector of matrix values
219 real*8 usrmem(*) ! optional user memory
220
221 Integer :: i, j, k, l, nzc
222!
223! Information about Variables:
224! Default: Lower = -Inf, Curr = 0, and Upper = +inf.
225! Default: the status information in Vsta is not used.
226!
227 l = 0
228 do i = 1, dimx
229 do k = 1, dimk
230 l = l + 1
231 curr(l) = +0.9d0 ! Make the initial point infeasible to the sum(j,(x(j,k)) =L= 1;
232 lower(l) = -1.0d0
233 enddo
234 enddo
235!
236! Information about Constraints:
237! Default: Rhs = 0
238! Default: the status information in Esta and the function
239! value in FV are not used.
240! Default: Type: There is no default.
241! 0 = Equality,
242! 1 = Greater than or equal,
243! 2 = Less than or equal,
244! 3 = Non binding.
245!
246! Constraints 1 to Nobs*Dimk:
247! Rhs = Obs(i,k) and type Equality
248!
249 l = 0
250 do i = 1, nobs
251 do k = 1, dimk
252 l = l + 1
253 rhs(l) = obs(i,k)
254 type(l) = 0
255 enddo
256 enddo
257!
258! Constraint Nobs*Dimk+1 to 2*Nobs*Dimk
259! Rhs = 0 (default) and type Less than or equal
260!
261 do i = nobs*dimk+1, 2*nobs*dimk
262 type(i) = 2
263 enddo
264!
265! Constraint 2*Nobs*Dimk+1 to 3*Nobs*Dimk
266! Rhs = 0 (default) and type Greater than or equal
267!
268 do i = 2*nobs*dimk+1, 3*nobs*dimk
269 type(i) = 1
270 enddo
271!
272! Constraint 3*Nobs*Dimk+1 to 3*Nobs*Dimk+Dimk
273! Rhs = 1 and type Less than or equal
274!
275 do i = 3*nobs*dimk+1, 3*nobs*dimk+dimk
276 rhs(i) = 1.0d0
277 type(i) = 2
278 enddo
279!
280! Objective, type Nonbinding
281!
282 type(3*nobs*dimk+dimk+1) = 3
283!
284! Information about the Jacobian. We use the standard method with
285! Rowno, Value, Nlflag and Colsta and we do not use Colno.
286!
287! Colsta = Start of column indices (No Defaults):
288! Rowno = Row indices
289! Value = Value of derivative (by default only linear
290! derivatives are used)
291! Nlflag = 0 for linear and 1 for nonlinear derivative
292! (not needed for completely linear models)
293!
294!
295! Indices
296! x(j,k) res(i,k) z(k)
297! i,k: NL L=1
298! i,k: L=1 nl and L=-1
299! i,k: L=1 nl and L=+1
300! k L=1
301! obj L=+1
302!
303! We map (i,k) -> k+DimK*(i-1)
304! and (j,k) -> k+DimK*(j-1)
305!
306 nzc = 1
307 do j = 1, dimx ! x(j,k)
308 do k = 1, dimk
309 l = k+dimk*(j-1)
310 colsta(l) = nzc
311 do i = 1, nobs
312 rowno(nzc) = k+dimk*(i-1)
313 nlflag(nzc) = 1
314 nzc = nzc + 1
315 enddo
316 rowno(nzc) = 3*nobs*dimk+k
317 nlflag(nzc) = 0
318 value(nzc) = 1.d0
319 nzc = nzc + 1
320 enddo
321 enddo
322 do i = 1, nobs ! res(i,k)
323 do k = 1, dimk
324 colsta(dimx*dimk+k+dimk*(i-1)) = nzc
325 rowno(nzc) = k+dimk*(i-1)
326 nlflag(nzc) = 0
327 value(nzc) = 1.d0
328 nzc = nzc + 1
329 rowno(nzc) = nobs*dimk+k+dimk*(i-1)
330 nlflag(nzc) = 0
331 value(nzc) = 1.d0
332 nzc = nzc + 1
333 rowno(nzc) = 2*nobs*dimk+k+dimk*(i-1)
334 nlflag(nzc) = 0
335 value(nzc) = 1.d0
336 nzc = nzc + 1
337 enddo
338 enddo
339 do k = 1, dimk ! z(k)
340 colsta(dimx*dimk+nobs*dimk+k) = nzc
341 do i = 1, nobs
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 if ( k < dimk ) then ! also nonlinear terms for z(k-1)
352 do i = 1, nobs
353 rowno(nzc) = nobs*dimk+k+1+dimk*(i-1)
354 nlflag(nzc) = 1
355 nzc = nzc + 1
356 rowno(nzc) = 2*nobs*dimk+k+1+dimk*(i-1)
357 nlflag(nzc) = 1
358 nzc = nzc + 1
359 enddo
360 endif
361 rowno(nzc) = 3*nobs*dimk+dimk+1 ! Objective
362 nlflag(nzc) = 0
363 value(nzc) = +1.d0
364 nzc = nzc + 1
365 enddo
366 colsta(dimx*dimk+nobs*dimk+dimk+1) = nzc
367
368 mm_readmatrix = 0 ! Return value means OK
369
370end Function mm_readmatrix
371!
372!==========================================================================
373! Compute nonlinear terms and non-constant Jacobian elements
374!
375
376!> Compute nonlinear terms and non-constant Jacobian elements
377!!
378!! @include{doc} fdeval_params.dox
379Integer Function mm_fdeval( x, g, jac, rowno, jcnm, mode, ignerr, errcnt, &
380 n, nzc, thread, usrmem )
381#if defined(itl)
382!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: MM_FDEval
383#endif
384 use lsq_7k
385 implicit none
386 integer, intent (in) :: n ! number of variables
387 integer, intent (in) :: rowno ! number of the row to be evaluated
388 integer, intent (in) :: nzc ! number of nonzceros in this row
389 real*8, intent (in), dimension(n) :: x ! vector of current solution values
390 real*8, intent (in out) :: g ! constraint value
391 real*8, intent (in out), dimension(n) :: jac ! vector of derivatives for current constraint
392 integer, intent (in), dimension(nzc) :: jcnm ! list of variables that appear nonlinearly
393 ! in this row. Ffor information only.
394 integer, intent (in) :: mode ! evaluation mode: 1 = function value
395 ! 2 = derivatives, 3 = both
396 integer, intent (in) :: ignerr ! if 1 then errors can be ignored as long
397 ! as errcnt is incremented
398 integer, intent (in out) :: errcnt ! error counter to be incremented in case
399 ! of function evaluation errors.
400 integer, intent (in) :: thread
401 real*8 usrmem(*) ! optional user memory
402
403 integer :: i,j,k,l
404 real*8 :: s
405
406 mm_fdeval = 0
407 if ( rowno .le. nobs*dimk ) then
408!
409! We map (i,k): rowno -> k+DimK*(i-1) so
410! i = (Dimk+rowno-1))/Dimk
411! k = rowno-Dimk*(i-1)
412! and (j,k) -> k+DimK*(j-1)
413! Mode = 1 or 3: Function value - x-part only
414!
415 i = (dimk+rowno-1)/dimk
416 k = rowno-dimk*(i-1)
417 if ( mode .eq. 1 .or. mode .eq. 3 ) then
418 s = 0.d0
419 do j = 1, dimx
420 l = k+dimk*(j-1)
421 s = s + a(i,k,j)*x(l) + b(i,k,j)*x(l)**2
422 enddo
423 g = s
424 endif
425!
426! Mode = 2 or 3: Derivatives
427!
428 if ( mode .eq. 2 .or. mode .eq. 3 ) then
429 do j = 1, dimx
430 l = k+dimk*(j-1)
431 jac(l) = a(i,k,j) + 2.d0*b(i,k,j)*x(l)
432 enddo
433 endif
434 Else if ( rowno .le. 2*nobs*dimk ) then
435 i = (dimk+rowno-nobs*dimk-1)/dimk
436 k = rowno-nobs*dimk-dimk*(i-1)-1
437 if ( mode .eq. 1 .or. mode .eq. 3 ) then
438 g = x(nobs*dimk + dimx*dimk+k )**2
439 endif
440 if ( mode .eq. 2 .or. mode .eq. 3 ) then
441 jac(nobs*dimk + dimx*dimk+k) = 2.0d0*x(nobs*dimk + dimx*dimk+k )
442 endif
443 Else if ( rowno .le. 3*nobs*dimk ) then
444 i = (dimk+rowno-2*nobs*dimk-1)/dimk
445 k = rowno-2*nobs*dimk-dimk*(i-1)-1
446 if ( mode .eq. 1 .or. mode .eq. 3 ) then
447 g = x(nobs*dimk + dimx*dimk+k )**2
448 endif
449 if ( mode .eq. 2 .or. mode .eq. 3 ) then
450 jac(nobs*dimk + dimx*dimk+k) = 2.0d0*x(nobs*dimk + dimx*dimk+k )
451 endif
452 Else
453 mm_fdeval = 1 ! Illegal row number
454 endif
455
456end 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 minimax07
Main program. A simple setup and call of CONOPT.
Definition minimax07.f90:85
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