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