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