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