CONOPT
Loading...
Searching...
No Matches
stress1.f90
Go to the documentation of this file.
1!> @file stress1.f90
2!! @ingroup FORT1THREAD_EXAMPLES
3!!
4!!
5!! This is a CONOPT implementation of the gams model
6!!
7!! @verbatim
8!! set i / i1*i10000/
9!! variable x(i), obj;
10!! equation objdef; objdef .. obj =E= prod(i, x(i)*x(i) );
11!! x.lo(i) = 1;
12!! model m / all /;
13!! solve m using nlp minimizing obj;
14!! @endverbatim
15!!
16!!
17!! For more information about the individual callbacks, please have a look at the source code.
18
19#if defined(_WIN32) && !defined(_WIN64)
20#define dec_directives_win32
21#endif
22
23Module stress1data
24 integer, Parameter :: N = 15000
25End Module stress1data
27!> Main program. A simple setup and call of CONOPT
28!!
29Program stress1
30 Use proginfo
32 Use stress1data
33 Implicit none
34!
35! Declare the user callback routines as Integer, External:
36!
37 Integer, External :: stress_readmatrix ! Mandatory Matrix definition routine defined below
38 Integer, External :: stress_fdeval ! Function and Derivative evaluation routine
39 ! needed a nonlinear model.
40 Integer, External :: std_status ! Standard callback for displaying solution status
41 Integer, External :: std_solution ! standard callback for displaying solution values
42 Integer, External :: std_message ! Standard callback for managing messages
43 Integer, External :: std_errmsg ! Standard callback for managing error messages
44 Integer, External :: stress_2dlagrstr ! Hessian structure evaluation routine
45 Integer, External :: stress_2dlagrval ! Hessian value evaluation routine
46#ifdef dec_directives_win32
47!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Stress_ReadMatrix
48!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Stress_FDEval
49!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Status
50!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Solution
51!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Message
52!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_ErrMsg
53!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Stress_2DLagrStr
54!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Stress_2DLagrVal
55#endif
56!
57! Control vector
58!
59 INTEGER, Dimension(:), Pointer :: cntvect
60 INTEGER :: coi_error
61
62 call startup
63!
64! Create and initialize a Control Vector
65!
66 coi_error = coi_create( cntvect )
67!
68! Number of variables (excl. slacks): N
69!
70 coi_error = max( coi_error, coidef_numvar( cntvect, n ) )
71!
72! Number of equations: 1 objective
73!
74 coi_error = max( coi_error, coidef_numcon( cntvect, 1 ) )
75!
76! Number of nonzeros in the Jacobian: N
77!
78 coi_error = max( coi_error, coidef_numnz( cntvect, n ) )
79!
80! Number of nonlinear nonzeros. N
81!
82 coi_error = max( coi_error, coidef_numnlnz( cntvect, n ) )
83!
84! Number of nonzeros in the Hessian. N*(N+1)/2
85!
86 coi_error = max( coi_error, coidef_numhess( cntvect, n*(n+1)/2 ) )
87!
88! Direction: -1 = minimization.
89!
90 coi_error = max( coi_error, coidef_optdir( cntvect, -1 ) )
91!
92! Objective: Constraint no 1
93!
94 coi_error = max( coi_error, coidef_objcon( cntvect, 1 ) )
95!
96! Turn Hessian debugging off
97!
98 coi_error = max( coi_error, coidef_debug2d( cntvect, 0 ) )
99!
100! Options file = stress1.opt
101!
102 coi_error = max( coi_error, coidef_optfile( cntvect, 'stress1.opt' ) )
103!
104! Tell CONOPT about the callback routines:
105!
106 coi_error = max( coi_error, coidef_readmatrix( cntvect, stress_readmatrix ) )
107 coi_error = max( coi_error, coidef_fdeval( cntvect, stress_fdeval ) )
108 coi_error = max( coi_error, coidef_status( cntvect, std_status ) )
109 coi_error = max( coi_error, coidef_solution( cntvect, std_solution ) )
110 coi_error = max( coi_error, coidef_message( cntvect, std_message ) )
111 coi_error = max( coi_error, coidef_errmsg( cntvect, std_errmsg ) )
112 coi_error = max( coi_error, coidef_2dlagrstr( cntvect, stress_2dlagrstr ) )
113 coi_error = max( coi_error, coidef_2dlagrval( cntvect, stress_2dlagrval ) )
114
115#if defined(CONOPT_LICENSE_INT_1) && defined(CONOPT_LICENSE_INT_2) && defined(CONOPT_LICENSE_INT_3) && defined(CONOPT_LICENSE_TEXT)
116 coi_error = max( coi_error, coidef_license( cntvect, conopt_license_int_1, conopt_license_int_2, conopt_license_int_3, conopt_license_text) )
117#endif
118
119 If ( coi_error .ne. 0 ) THEN
120 write(*,*)
121 write(*,*) '**** Fatal Error while loading CONOPT Callback routines.'
122 write(*,*)
123 call flog( "Skipping First Solve due to setup errors", 1 )
124 ENDIF
125!
126! Start CONOPT:
127!
128 coi_error = coi_solve( cntvect )
129 If ( coi_error .ne. 0 ) THEN
130 write(*,*)
131 write(*,*) '**** Fatal Error while Solving model. Return Code=',coi_error
132 write(*,*)
133 call flog( "Errors encountered during Solve.", 1 )
134 elseif ( stacalls == 0 .or. solcalls == 0 ) then
135 call flog( "Status or Solution routine was not called during Solve", 1 )
136 else if ( mstat /= 2 .or. sstat /= 1 ) then
137 call flog( "Incorrect Model or Solver Status during solve", 1 )
138 ENDIF
139
140 write(*,*)
141 write(*,*) 'End of Stress1 Model. Return code=',coi_error
142
143 if ( coi_free(cntvect) /= 0 ) call flog( "Error while freeing control vector",1)
144
145 call flog( "Successful Solve", 0 )
146!
147! Free solution memory
148!
150
151end Program stress1
152!
153! =====================================================================
154! Define information about the model structure
155!
156
157!> Define information about the model
158!!
159!! @include{doc} readMatrix_params.dox
160Integer Function stress_readmatrix( lower, curr, upper, vsta, type, rhs, esta, &
161 colsta, rowno, value, nlflag, nc, mc, nz, usrmem )
162#ifdef dec_directives_win32
163!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Stress_ReadMatrix
164#endif
165 Use stress1data
166 implicit none
167 integer, intent (in) :: nc ! number of variables
168 integer, intent (in) :: mc ! number of constraints
169 integer, intent (in) :: nz ! number of nonzeros
170 real*8, intent (in out), dimension(nc) :: lower ! vector of lower bounds
171 real*8, intent (in out), dimension(nc) :: curr ! vector of initial values
172 real*8, intent (in out), dimension(nc) :: upper ! vector of upper bounds
173 integer, intent (in out), dimension(nc) :: vsta ! vector of initial variable status
174 ! Defined during restarts
175 integer, intent (out), dimension(mc) :: type ! vector of equation types
176 integer, intent (in out), dimension(mc) :: esta ! vector of initial equation status
177 ! Defined during restarts
178 real*8, intent (in out), dimension(mc) :: rhs ! vector of right hand sides
179 integer, intent (in out), dimension(nc+1) :: colsta ! vector with start of column indices
180 integer, intent (out), dimension(nz) :: rowno ! vector of row numbers
181 integer, intent (in out), dimension(nz) :: nlflag ! vector of nonlinearity flags
182 real*8, intent (in out), dimension(nz) :: value ! vector of matrix values
183 real*8 usrmem(*) ! optional user memory
184
185 Integer :: i
186!
187! Define the information for the columns.
188!
189! We should not supply status information, vsta.
190!
191 do i = 1, n
192 lower(i) = 1.d0
193 curr(i) = 1.d0
194 enddo
195!
196! Define the information for the rows
197!
198! We order the constraints as follows: The objective is first,
199! followed by tddef, sdef, csdef, ddef, rdef, and revdef for
200! the first period, the same for the second period, etc.
201!
202! The objective is a nonbinding constraint:
203!
204 type(1) = 3
205!
206! Define the structure and content of the Jacobian:
207!
208 do i = 1, n
209 colsta(i) = i
210 rowno(i) = 1
211 nlflag(i) = 1
212 enddo
213 colsta(n+1) = n+1
214
216
217end Function stress_readmatrix
218!
219! =====================================================================
220! Compute nonlinear terms and non-constant Jacobian elements
221!
222
223!> Compute nonlinear terms and non-constant Jacobian elements
224!!
225!! @include{doc} fdeval_params.dox
226Integer Function stress_fdeval( x, g, jac, rowno, jcnm, mode, ignerr, errcnt, &
227 nc, nz, thread, usrmem )
228#ifdef dec_directives_win32
229!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Stress_FDEval
230#endif
231 Use stress1data
232 implicit none
233 integer, intent (in) :: nc ! number of variables
234 integer, intent (in) :: rowno ! number of the row to be evaluated
235 integer, intent (in) :: nz ! number of nonzeros in this row
236 real*8, intent (in), dimension(nc) :: x ! vector of current solution values
237 real*8, intent (in out) :: g ! constraint value
238 real*8, intent (in out), dimension(nc) :: jac ! vector of derivatives for current constraint
239 integer, intent (in), dimension(nz) :: jcnm ! list of variables that appear nonlinearly
240 ! in this row. Ffor information only.
241 integer, intent (in) :: mode ! evaluation mode: 1 = function value
242 ! 2 = derivatives, 3 = both
243 integer, intent (in) :: ignerr ! if 1 then errors can be ignored as long
244 ! as errcnt is incremented
245 integer, intent (in out) :: errcnt ! error counter to be incremented in case
246 ! of function evaluation errors.
247 integer, intent (in) :: thread
248 real*8 usrmem(*) ! optional user memory
249
250 integer :: i
251 real*8 :: obj
252!
253 if ( rowno == 1 ) then
254 obj = 1.0d0
255 do i = 1, n
256 obj = obj * x(i)*x(i)
257 enddo
258 if ( mode == 1 .or. mode == 3 ) then
259 g = obj
260 endif
261 if ( mode == 2 .or. mode == 3 ) then
262 do i = 1, n
263 jac(i) = 2.d0*obj/x(i)
264 enddo
265 endif
266 else
268! Error - this equation is not nonlinear
269!
270 write(*,*)
271 write(*,*) 'Error. FDEval called with rowno=',rowno
272 write(*,*)
273 stress_fdeval = 1
274 Return
275 endif
276 stress_fdeval = 0
277
278end Function stress_fdeval
279
280
281!> Specify the structure of the Lagrangian of the Hessian
282!!
283!! @include{doc} 2DLagrStr_params.dox
284Integer Function stress_2dlagrstr( HSRW, HSCL, NODRV, Nc, Mc, NHESS, UsrMem )
285#ifdef dec_directives_win32
286!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Stress_2DLagrStr
287#endif
288 USE stress1data
289 Implicit None
290
291 Integer, Intent (IN) :: nc, mc, nhess
292 Integer, Intent (IN OUT) :: nodrv
293 Integer, Dimension(NHess), Intent (Out) :: hsrw, hscl
294 real*8, Intent(IN OUT) :: usrmem(*)
295
296 Integer :: i, j, iz
297!
298! Define structure of Hessian
299!
300 if ( n /= nc ) then
301 write(*,*) 'Struct: Number of variables should be=',n,' but CONOPT uses',nc
302 stress_2dlagrstr = 1; Return
303 stop
304 endif
305 if ( mc /= 1 ) then
306 write(*,*) 'Struct: Number of equations should be=',+1,' but CONOPT uses',mc
307 stress_2dlagrstr = 1; Return
308 endif
309 if ( nhess /= n*(n+1)/2 ) then
310 write(*,*) 'Struct: Number of Hessian elements should be=',n*(n+1)/2,' but CONOPT uses',nhess
311 stress_2dlagrstr = 1; Return
312 endif
313 iz = 0
314 Do i = 1, n
315 Do j = i, n
316 iz = iz + 1
317 hscl(iz) = i
318 hsrw(iz) = j
319 Enddo
320 Enddo
321
323
324End Function stress_2dlagrstr
325
326
327!> Compute the Lagrangian of the Hessian
328!!
329!! @include{doc} 2DLagrVal_params.dox
330Integer Function stress_2dlagrval( X, U, HSRW, HSCL, HSVL, NODRV, Nc, Mc, NHESS, UsrMem )
331#ifdef dec_directives_win32
332!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Stress_2DLagrVal
333#endif
334 USE stress1data
335 Implicit None
336
337 Integer, Intent (IN) :: nc, mc, nhess
338 Integer, Intent (IN OUT) :: nodrv
339 real*8, Dimension(Nc), Intent (IN) :: x
340 real*8, Dimension(Mc), Intent (IN) :: u
341 Integer, Dimension(Nhess), Intent (IN) :: hsrw, hscl
342 real*8, Dimension(NHess), Intent (Out) :: hsvl
343 real*8, Intent(IN OUT) :: usrmem(*)
344
345 Integer :: i, j, iz
346!
347! Normal Evaluation mode -- we are testing memory only so we return a unit matrix.
348!
349 iz = 0
350 Do i = 1, n
351 iz = iz + 1
352 hsvl(iz) = 1.0d0
353 Do j = i+1, n
354 iz = iz + 1
355 hsvl(iz) = 0.0d0
356 Enddo
357 Enddo
358
360
361End Function stress_2dlagrval
integer function std_solution(xval, xmar, xbas, xsta, yval, ymar, ybas, ysta, n, m, usrmem)
Definition comdecl.f90:170
integer function std_status(modsta, solsta, iter, objval, usrmem)
Definition comdecl.f90:126
integer function std_message(smsg, dmsg, nmsg, llen, usrmem, msgv)
Definition comdecl.f90:243
integer function std_errmsg(rowno, colno, posno, msglen, usrmem, msg)
Definition comdecl.f90:286
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_2dlagrstr(cntvect, coi_2dlagrstr)
define callback routine for providing the structure of the second derivatives of the Lagrangian.
Definition conopt.f90:1547
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_2dlagrval(cntvect, coi_2dlagrval)
define callback routine for computing the values of the second derivatives of the Lagrangian.
Definition conopt.f90:1573
integer(c_int) function coidef_license(cntvect, licint1, licint2, licint3, licstring)
define the License Information.
Definition conopt.f90:293
integer(c_int) function coidef_debug2d(cntvect, debug2d)
turn debugging of 2nd derivatives on and off.
Definition conopt.f90:534
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_numhess(cntvect, numhess)
defines the Number of Hessian Nonzeros.
Definition conopt.f90:188
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
integer solcalls
Definition comdecl.f90:15
integer sstat
Definition comdecl.f90:18
subroutine finalize
Definition comdecl.f90:79
integer stacalls
Definition comdecl.f90:14
subroutine flog(msg, code)
Definition comdecl.f90:62
integer mstat
Definition comdecl.f90:17
subroutine startup
Definition comdecl.f90:41
integer, parameter n
Definition stress1.f90:26
integer function stress_2dlagrstr(hsrw, hscl, nodrv, nc, mc, nhess, usrmem)
Specify the structure of the Lagrangian of the Hessian.
Definition stress1.f90:268
integer function stress_readmatrix(lower, curr, upper, vsta, type, rhs, esta, colsta, rowno, value, nlflag, nc, mc, nz, usrmem)
Define information about the model.
Definition stress1.f90:151
integer function stress_2dlagrval(x, u, hsrw, hscl, hsvl, nodrv, nc, mc, nhess, usrmem)
Compute the Lagrangian of the Hessian.
Definition stress1.f90:311
program stress1
Main program. A simple setup and call of CONOPT.
Definition stress1.f90:31
integer function stress_fdeval(x, g, jac, rowno, jcnm, mode, ignerr, errcnt, nc, nz, thread, usrmem)
Compute nonlinear terms and non-constant Jacobian elements.
Definition stress1.f90:214