CONOPT
Loading...
Searching...
No Matches
mp_rosex.f90
Go to the documentation of this file.
1!> @file mp_rosex.f90
2!! @ingroup FORTOPENMP_EXAMPLES
3!!
4!! @copydoc rosex.f90
5
6Module tlim
7 Integer :: thread2d
8End module tlim
9
10!> Main program. A simple setup and call of CONOPT
11!!
12Program rosex
13
14 Use proginfop
15 Use coidef
16 Use omp_lib
17 Use tlim
18 IMPLICIT NONE
19!
20! Declare the user callback routines as Integer, External:
21!
22 Integer, External :: ros_readmatrix ! Mandatory Matrix definition routine defined below
23 Integer, External :: ros_fdeval ! Function and Derivative evaluation routine
24 ! needed a nonlinear model.
25 Integer, External :: std_status ! Standard callback for displaying solution status
26 Integer, External :: ros_solution ! Special callback for displaying solution values that does not write the solution
27 Integer, External :: std_message ! Standard callback for managing messages
28 Integer, External :: std_errmsg ! Standard callback for managing error messages
29#if defined(itl)
30!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Ros_ReadMatrix
31!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Ros_FDEval
32!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Status
33!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Ros_Solution
34!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Message
35!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_ErrMsg
36#endif
37!
38! Control vector
39!
40 INTEGER, Dimension(:), Pointer :: cntvect
41 INTEGER :: coi_error
42!
43! The model is an unconstrained problem in which the objective is the
44! sum of N identical copies of the Rosenbrock function placed after
45! each other:
46!
47! min sum(i, (1-x(2*i-1))**2+100*(x(2*i)-x(2*i-1)**2)**2)
48!
49! The number of Rosenborck functions
50!
51 Integer :: nr
52!
53! Thread Info
54!
55 Integer :: maxthread
56 real*8 time0, time1, time4
57!
58! Create and initialize a Control Vector
59!
60 call startup
61
62 coi_error = coi_createfort( cntvect )
63!
64! Tell CONOPT about the size of the model by populating the Control Vector:
65!
66 nr = 16000
67 coi_error = max( coi_error, coidef_numvar( cntvect, 2*nr ) ) ! # variables
68 coi_error = max( coi_error, coidef_numcon( cntvect, 1 ) ) ! # constraints
69 coi_error = max( coi_error, coidef_numnz( cntvect, 2*nr ) ) ! # nonzeros in the Jacobian
70 coi_error = max( coi_error, coidef_numnlnz( cntvect, 2*nr ) ) ! # of which are nonlinear
71 coi_error = max( coi_error, coidef_optdir( cntvect, -1 ) ) ! Minimize
72 coi_error = max( coi_error, coidef_objcon( cntvect, 1 ) ) ! Objective is constraint 1
73 coi_error = max( coi_error, coidef_optfile( cntvect, 'mp_rosex.opt' ) )
74!
75! Tell CONOPT about the callback routines:
76!
77 coi_error = max( coi_error, coidef_readmatrix( cntvect, ros_readmatrix ) )
78 coi_error = max( coi_error, coidef_fdeval( cntvect, ros_fdeval ) )
79 coi_error = max( coi_error, coidef_status( cntvect, std_status ) )
80 coi_error = max( coi_error, coidef_solution( cntvect, ros_solution ) )
81 coi_error = max( coi_error, coidef_message( cntvect, std_message ) )
82 coi_error = max( coi_error, coidef_errmsg( cntvect, std_errmsg ) )
83
84#if defined(LICENSE_INT_1) && defined(LICENSE_INT_2) && defined(LICENSE_INT_3) && defined(LICENSE_TEXT)
85 coi_error = max( coi_error, coidef_license( cntvect, license_int_1, license_int_2, license_int_3, license_text) )
86#endif
87
88 If ( coi_error .ne. 0 ) THEN
89 write(*,*)
90 write(*,*) '**** Fatal Error while loading CONOPT Callback routines.'
91 write(*,*)
92 call flog( "Skipping Solve due to setup errors", 1 )
93 ENDIF
94!
95! Start CONOPT with a single thread:
96!
97 thread2d = 1
98 time0 = omp_get_wtime()
99 coi_error = coi_solve( cntvect )
100 time1 = omp_get_wtime() - time0
101 If ( coi_error /= 0 ) then
102 call flog( "Solve 1: Errors encountered during solution", 1 )
103 elseif ( stacalls == 0 .or. solcalls == 0 ) then
104 call flog( "Solve 1: Status or Solution routine was not called", 1 )
105 elseif ( sstat /= 1 .or. mstat /= 2 ) then
106 call flog( "Solve 1: Solver and Model Status was not as expected (1,2)", 1 )
107 elseif ( abs( obj-0.0d0 ) > 0.000001d0 ) then
108 call flog( "Solve 1: Incorrect objective returned", 1 )
109 endif
110!
111 maxthread = omp_get_max_threads() ! Define ThreadC as the largest number of threads we can get
112 coi_error = max( coi_error, coidef_threadc( cntvect, maxthread ) )
113 coi_error = max( coi_error, coidef_threadf( cntvect, 1 ) ) ! Functions must run with one thread
114#if defined (complete)
115!
116! Start CONOPT -- Multi-thread mode with a single thread:
117!
118 thread2d = 1
119 coi_error = max( coi_error, coidef_threads( cntvect, 1 ) ) ! 1 means use 1 thread
120 time0 = omp_get_wtime()
121 coi_error = coi_solve( cntvect )
122 time1a = omp_get_wtime() - time0
123 If ( coi_error /= 0 ) then
124 call flog( "Solve 1a: Errors encountered during solution", 1 )
125 elseif ( stacalls == 0 .or. solcalls == 0 ) then
126 call flog( "Solve 1a: Status or Solution routine was not called", 1 )
127 elseif ( sstat /= 1 .or. mstat /= 2 ) then
128 call flog( "Solve 1a: Solver and Model Status was not as expected (1,2)", 1 )
129 elseif ( abs( obj-0.0d0 ) > 0.000001d0 ) then
130 call flog( "Solve 1a: Incorrect objective returned", 1 )
131 endif
132!
133! Start CONOPT -- Multi-thread mode with two threads:
134!
135 thread2d = 2
136 coi_error = max( coi_error, coidef_threads( cntvect, 2 ) ) ! 2 threads
137 time0 = omp_get_wtime()
138 coi_error = coi_solve( cntvect )
139 time2 = omp_get_wtime() - time0
140 If ( coi_error /= 0 ) then
141 call flog( "Solve 2: Errors encountered during solution", 1 )
142 elseif ( stacalls == 0 .or. solcalls == 0 ) then
143 call flog( "Solve 2: Status or Solution routine was not called", 1 )
144 elseif ( sstat /= 1 .or. mstat /= 2 ) then
145 call flog( "Solve 2: Solver and Model Status was not as expected (1,2)", 1 )
146 elseif ( abs( obj-0.0d0 ) > 0.000001d0 ) then
147 call flog( "Solve 2: Incorrect objective returned", 1 )
148 endif
149#endif
150 if ( maxthread >= 4 ) then
151!
152! Start CONOPT -- Multi-thread mode with four threads:
153!
154 thread2d = 4
155 coi_error = max( coi_error, coidef_threads( cntvect, 4 ) ) ! 4 threads
156 time0 = omp_get_wtime()
157 coi_error = coi_solve( cntvect )
158 time4 = omp_get_wtime() - time0
159 If ( coi_error /= 0 ) then
160 call flog( "Solve 4: Errors encountered during solution", 1 )
161 elseif ( stacalls == 0 .or. solcalls == 0 ) then
162 call flog( "Solve 4: Status or Solution routine was not called", 1 )
163 elseif ( sstat /= 1 .or. mstat /= 2 ) then
164 call flog( "Solve 4: Solver and Model Status was not as expected (1,2)", 1 )
165 elseif ( abs( obj-0.0d0 ) > 0.000001d0 ) then
166 call flog( "Solve 4: Incorrect objective returned", 1 )
167 endif
168 endif
169
170 if ( coi_free( cntvect ) /= 0 ) call flog( "Error while freeing control vector", 1 )
171
172 write(*,*)
173 write(*,"('Time for single thread ',f10.3)") time1
174#if defined (complete)
175 write(*,"('Time for single thread ',f10.3)") time1a
176 write(*,"('Time for dual threads ',f10.3)") time2
177#endif
178 if ( maxthread >= 4 ) then
179 write(*,"('Time for quad threads ',f10.3)") time4
180 endif
181 write(*,*)
182 write(*,*) 'End of Extended Rosenbrock Function example. Return code=',coi_error
183
184 call flog( "Successful Solve", 0 )
185
186End Program rosex
187!
188! ============================================================================
189! Define information about the model:
190!
191
192!> Define information about the model
193!!
194!! @include{doc} readMatrix_params.dox
195Integer Function ros_readmatrix( lower, curr, upper, vsta, type, rhs, esta, &
196 colsta, rowno, value, nlflag, n, m, nz, &
197 usrmem )
198#if defined(itl)
199!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Ros_ReadMatrix
200#endif
201 IMPLICIT NONE
202 integer, intent (in) :: n ! number of variables
203 integer, intent (in) :: m ! number of constraints
204 integer, intent (in) :: nz ! number of nonzeros
205 real*8, intent (in out), dimension(n) :: lower ! vector of lower bounds
206 real*8, intent (in out), dimension(n) :: curr ! vector of initial values
207 real*8, intent (in out), dimension(n) :: upper ! vector of upper bounds
208 integer, intent (in out), dimension(n) :: vsta ! vector of initial variable status
209 ! (not defined here)
210 integer, intent (out), dimension(m) :: type ! vector of equation types
211 integer, intent (in out), dimension(m) :: esta ! vector of initial equation status
212 ! (not defined here)
213 real*8, intent (in out), dimension(m) :: rhs ! vector of right hand sides
214 integer, intent (in out), dimension(n+1) :: colsta ! vector with start of column indices
215 integer, intent (out), dimension(nz) :: rowno ! vector of row numbers
216 integer, intent (in out), dimension(nz) :: nlflag ! vector of nonlinearity flags
217 real*8, intent (in out), dimension(nz) :: value ! vector of matrix values
218 real*8 usrmem(*) ! optional user memory
219
220 Integer :: nr
221 Integer :: i
222
223 nr = n / 2
224!
225! Information about Variables:
226! Default: Lower = -Inf, Curr = 0, and Upper = +inf.
227! Default: the status information in Vsta is not used.
228!
229 DO i = 1, nr
230 curr(2*i-1) = -1.d0
231 curr(2*i ) = +1.d0
232 enddo
233!
234! Information about Constraints:
235! Default: Rhs = 0
236! Default: the status information in Esta and the function
237! value in FV are not used.
238! Default: Type: There is no default.
239! 0 = Equality,
240! 1 = Greater than or equal,
241! 2 = Less than or equal,
242! 3 = Non binding.
243!
244! Constraint 1 (Objective)
245! Rhs = 0.0 and type Non binding
246!
247 type(1) = 3
248!
249! Information about the Jacobian. We have to define Rowno, Value,
250! Nlflag and Colsta.
251!
252! Colsta = Start of column indices (No Defaults):
253! Rowno = Row indices
254! Value = Value of derivative (by default only linear
255! derivatives are used)
256! Nlflag = 0 for linear and 1 for nonlinear derivative
257! (not needed for completely linear models)
258!
259! Indices
260! x(1) x(2) x(3) x(4)
261! 1: 1 2 3 4 etc
262!
263! Nonlinearity Structure: L = 0 are linear and NL = 1 are nonlinear
264! x(1) x(2) x(3) x(4)
265! 1: NL NL NL NL etc
266!
267! Value (Linear only)
268!
269 DO i = 1, nr*2+1
270 colsta(i) = i
271 enddo
272 DO i = 1, nr*2
273 rowno(i) = 1
274 nlflag(i) = 1
275 enddo
276
277 ros_readmatrix = 0 ! Return value means OK
278
279end Function ros_readmatrix
280!
281!==========================================================================
282! Compute nonlinear terms and non-constant Jacobian elements
283!
284
285!> Compute nonlinear terms and non-constant Jacobian elements
286!!
287!! @include{doc} fdeval_params.dox
288Integer Function ros_fdeval( x, g, jac, rowno, jcnm, mode, ignerr, errcnt, &
289 n, nz, thread, usrmem )
290#if defined(itl)
291!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Ros_FDEval
292#endif
293 Use tlim
294 IMPLICIT NONE
295 integer, intent (in) :: n ! number of variables
296 integer, intent (in) :: rowno ! number of the row to be evaluated
297 integer, intent (in) :: nz ! number of nonzeros in this row
298 real*8, intent (in), dimension(n) :: x ! vector of current solution values
299 real*8, intent (in out) :: g ! constraint value
300 real*8, intent (in out), dimension(n) :: jac ! vector of derivatives for current constraint
301 integer, intent (in), dimension(nz) :: jcnm ! list of variables that appear nonlinearly
302 ! in this row. Ffor information only.
303 integer, intent (in) :: mode ! evaluation mode: 1 = function value
304 ! 2 = derivatives, 3 = both
305 integer, intent (in) :: ignerr ! if 1 then errors can be ignored as long
306 ! as errcnt is incremented
307 integer, intent (in out) :: errcnt ! error counter to be incremented in case
308 ! of function evaluation errors.
309 integer, intent (in) :: thread
310 real*8 usrmem(*) ! optional user memory
311
312 Integer :: nr
313 Integer :: i
314
315 nr = n / 2
316!
317! Row 1: the objective function is nonlinear
318!
319 if ( rowno .eq. 1 ) then
320!
321! Mode = 1 or 3. Function value:
322!
323 if ( mode .eq. 1 .or. mode .eq. 3 ) then
324 g = 0.d0
325 do i = 1, nr
326 g = g + (x(2*i-1)-1.d0)**2 + 1.d2*(x(2*i)-x(2*i-1)**2)**2
327 enddo
328 endif
329!
330! Mode = 2 or 3: Derivative values:
331!
332 if ( mode .eq. 2 .or. mode .eq. 3 ) then
333!$OMP PARALLEL DEFAULT(Shared) IF( Thread2D > 1 ) NUM_THREADS(Thread2D)
334!$OMP DO SCHEDULE(Static)
335 do i = 1, nr
336 jac(2*i-1) = 2.d0*(x(2*i-1)-1.d0) - 4.d2*(x(2*i)-x(2*i-1)**2)*x(2*i-1)
337 jac(2*i ) = 2.d2*(x(2*i)-x(2*i-1)**2)
338 enddo
339!$OMP END DO
340!$OMP END PARALLEL
341 endif
342 endif
343 ros_fdeval = 0
344
345end Function ros_fdeval
346
347Integer Function ros_solution( XVAL, XMAR, XBAS, XSTA, YVAL, YMAR, YBAS, YSTA, N, M, USRMEM )
348#if defined(itl)
349!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Ros_Solution
350#endif
351!
352! Simple implementation in which we write the solution values to
353! the 'Documentation file' on unit 10.
354!
355 Use proginfop
356 IMPLICIT NONE
357 INTEGER, Intent(IN) :: n, m
358 INTEGER, Intent(IN), Dimension(N) :: xbas, xsta
359 INTEGER, Intent(IN), Dimension(M) :: ybas, ysta
360 real*8, Intent(IN), Dimension(N) :: xval, xmar
361 real*8, Intent(IN), Dimension(M) :: yval, ymar
362 real*8, Intent(IN OUT) :: usrmem(*)
363
364
365 solcalls = solcalls + 1
366 ros_solution = 0
367
368END Function ros_solution
369
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_optfile(cntvect, optfile)
define callback routine for defining an options file.
integer function coidef_threadf(cntvect, threadf)
number of threads allowed for simultaneous FDEval calls.
integer function coidef_threads(cntvect, threads)
number of threads allowed internally in CONOPT.
integer function coidef_license(cntvect, licint1, licint2, licint3, licstring)
define the License Information.
Definition coistart.f90:680
integer function coidef_threadc(cntvect, threadc)
check for thread compatibility.
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 coi_solve(cntvect)
method for starting the solving process of CONOPT.
Definition coistart.f90:14
integer function ros_solution(xval, xmar, xbas, xsta, yval, ymar, ybas, ysta, n, m, usrmem)
Definition mp_rosex.f90:348
subroutine flog(msg, code)
Definition comdeclp.f90:42
real *8 obj
Definition comdeclp.f90:13
subroutine startup
Definition comdeclp.f90:25
integer thread2d
Definition mp_rosex.f90:7
integer function ros_readmatrix(lower, curr, upper, vsta, type, rhs, esta, colsta, rowno, value, nlflag, n, m, nz, usrmem)
Define information about the model.
Definition roseq.f90:120
integer function ros_fdeval(x, g, jac, rowno, jcnm, mode, ignerr, errcnt, n, nz, thread, usrmem)
Compute nonlinear terms and non-constant Jacobian elements.
Definition roseq.f90:211
program rosex
Main program. A simple setup and call of CONOPT.
Definition rosex.f90:12