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