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!! @verbatim
11!! min sum(k, max(i, abs(res(i,k)) )
12!!
13!! sum(j, a(i,j,k)*x(j,k) + b(i,j,k)*x(j,k)**2 ) + res(i,k) = obs(i,k)
14!! all k: sum(j,x(j,k)) =L= 1;
15!! @endverbatim
16!!
17!! where a, b, and obs are known data, and res and x are the
18!! variables of the model.
19!!
20!! To get a model with continuous derivatives we introduce an
21!! objective variable, z, and the constraints
22!! @verbatim
23!! +res(i,k) <= z(k) or +res(i,k) - z(k) <= 0
24!! -res(i,k) <= z(k) or +res(i,k) + z(k) >= 0
25!! @endverbatim
26!! Note opposite sign from Minimax01
27!! 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 minimax05
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, 'Minimax05.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 Minimax05 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( 'Minimax05', 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 minimax05
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.8d0
224 enddo
225 enddo
226!
227! Information about Constraints:
228! Default: Rhs = 0
229! Default: the status information in Esta and the function
230! value in FV are not used.
231! Default: Type: There is no default.
232! 0 = Equality,
233! 1 = Greater than or equal,
234! 2 = Less than or equal,
235! 3 = Non binding.
236!
237! Constraints 1 to Nobs*Dimk:
238! Rhs = Obs(i,k) and type Equality
239!
240 l = 0
241 do i = 1, nobs
242 do k = 1, dimk
243 l = l + 1
244 rhs(l) = obs(i,k)
245 type(l) = 0
246 enddo
247 enddo
248!
249! Constraint Nobs*Dimk+1 to 2*Nobs*Dimk
250! Rhs = 0 (default) and type Less than or equal
251!
252 do i = nobs*dimk+1, 2*nobs*dimk
253 type(i) = 2
254 enddo
255!
256! Constraint 2*Nobs*Dimk+1 to 3*Nobs*Dimk
257! Rhs = 0 (default) and type Greater than or equal
258!
259 do i = 2*nobs*dimk+1, 3*nobs*dimk
260 type(i) = 1
261 enddo
262!
263! Constraint 3*Nobs*Dimk+1 to 3*Nobs*Dimk+Dimk
264! Rhs = 1 and type Less than or equal
265!
266 do i = 3*nobs*dimk+1, 3*nobs*dimk+dimk
267 rhs(i) = 1.0d0
268 type(i) = 2
269 enddo
270!
271! Objective, type Nonbinding
272!
273 type(3*nobs*dimk+dimk+1) = 3
274!
275! Information about the Jacobian. We use the standard method with
276! Rowno, Value, Nlflag and Colsta and we do not use Colno.
277!
278! Colsta = Start of column indices (No Defaults):
279! Rowno = Row indices
280! Value = Value of derivative (by default only linear
281! derivatives are used)
282! Nlflag = 0 for linear and 1 for nonlinear derivative
283! (not needed for completely linear models)
284!
285!
286! Indices
287! x(j,k) res(i,k) z(k)
288! i,k: NL L=1
289! i,k: L=1 L=-1
290! i,k: L=1 L=+1
291! k L=1
292! obj L=+1
293!
294! We map (i,k) -> k+DimK*(i-1)
295! and (j,k) -> k+DimK*(j-1)
296!
297 nzc = 1
298 do j = 1, dimx ! x(j,k)
299 do k = 1, dimk
300 l = k+dimk*(j-1)
301 colsta(l) = nzc
302 do i = 1, nobs
303 rowno(nzc) = k+dimk*(i-1)
304 nlflag(nzc) = 1
305 nzc = nzc + 1
306 enddo
307 rowno(nzc) = 3*nobs*dimk+k
308 nlflag(nzc) = 0
309 value(nzc) = 1.d0
310 nzc = nzc + 1
311 enddo
312 enddo
313 do i = 1, nobs ! res(i,k)
314 do k = 1, dimk
315 colsta(dimx*dimk+k+dimk*(i-1)) = nzc
316 rowno(nzc) = k+dimk*(i-1)
317 nlflag(nzc) = 0
318 value(nzc) = 1.d0
319 nzc = nzc + 1
320 rowno(nzc) = nobs*dimk+k+dimk*(i-1)
321 nlflag(nzc) = 0
322 value(nzc) = 1.d0
323 nzc = nzc + 1
324 rowno(nzc) = 2*nobs*dimk+k+dimk*(i-1)
325 nlflag(nzc) = 0
326 value(nzc) = 1.d0
327 nzc = nzc + 1
328 enddo
329 enddo
330 do k = 1, dimk ! z(k)
331 colsta(dimx*dimk+nobs*dimk+k) = nzc
332 do i = 1, nobs
333 rowno(nzc) = nobs*dimk++k+dimk*(i-1)
334 nlflag(nzc) = 0
335 value(nzc) = -1.d0
336 nzc = nzc + 1
337 rowno(nzc) = 2*nobs*dimk+k+dimk*(i-1)
338 nlflag(nzc) = 0
339 value(nzc) = +1.d0
340 nzc = nzc + 1
341 enddo
342 rowno(nzc) = 3*nobs*dimk+dimk+1 ! Objective
343 nlflag(nzc) = 0
344 value(nzc) = +1.d0
345 nzc = nzc + 1
346 enddo
347 colsta(dimx*dimk+nobs*dimk+dimk+1) = nzc
348
349 mm_readmatrix = 0 ! Return value means OK
350
351end Function mm_readmatrix
352!
353!==========================================================================
354! Compute nonlinear terms and non-constant Jacobian elements
355!
356
357!> Compute nonlinear terms and non-constant Jacobian elements
358!!
359!! @include{doc} fdeval_params.dox
360Integer Function mm_fdeval( x, g, jac, rowno, jcnm, mode, ignerr, errcnt, &
361 n, nzc, thread, usrmem )
362#if defined(itl)
363!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: MM_FDEval
364#endif
365 use lsq_7k
366 implicit none
367 integer, intent (in) :: n ! number of variables
368 integer, intent (in) :: rowno ! number of the row to be evaluated
369 integer, intent (in) :: nzc ! number of nonzceros in this row
370 real*8, intent (in), dimension(n) :: x ! vector of current solution values
371 real*8, intent (in out) :: g ! constraint value
372 real*8, intent (in out), dimension(n) :: jac ! vector of derivatives for current constraint
373 integer, intent (in), dimension(nzc) :: jcnm ! list of variables that appear nonlinearly
374 ! in this row. Ffor information only.
375 integer, intent (in) :: mode ! evaluation mode: 1 = function value
376 ! 2 = derivatives, 3 = both
377 integer, intent (in) :: ignerr ! if 1 then errors can be ignored as long
378 ! as errcnt is incremented
379 integer, intent (in out) :: errcnt ! error counter to be incremented in case
380 ! of function evaluation errors.
381 integer, intent (in) :: thread
382 real*8 usrmem(*) ! optional user memory
383
384 integer :: i,j,k,l
385 real*8 :: s
386
387 mm_fdeval = 0
388 if ( rowno .le. nobs*dimk ) then
389!
390! We map (i,k): rowno -> k+DimK*(i-1) so
391! i = (Dimk+rowno-1))/Dimk
392! k = rowno-Dimk*(i-1)
393! and (j,k) -> k+DimK*(j-1)
394! Mode = 1 or 3: Function value - x-part only
395!
396 i = (dimk+rowno-1)/dimk
397 k = rowno-dimk*(i-1)
398 if ( mode .eq. 1 .or. mode .eq. 3 ) then
399 s = 0.d0
400 do j = 1, dimx
401 l = k+dimk*(j-1)
402 s = s + a(i,k,j)*x(l) + b(i,k,j)*x(l)**2
403 enddo
404 g = s
405 endif
406!
407! Mode = 2 or 3: Derivatives
408!
409 if ( mode .eq. 2 .or. mode .eq. 3 ) then
410 do j = 1, dimx
411 l = k+dimk*(j-1)
412 jac(l) = a(i,k,j) + 2.d0*b(i,k,j)*x(l)
413 enddo
414 endif
415 Else
416 mm_fdeval = 1 ! Illegal row number
417 endif
418
419end 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 minimax05
Main program. A simple setup and call of CONOPT.
Definition minimax05.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