CONOPT
Loading...
Searching...
No Matches
mp_lincns2.f90
Go to the documentation of this file.
1!> @file mp_lincns2.f90
2!! @ingroup FORTOPENMP_EXAMPLES
3!!
4!!
5!! Large Linear Dense CNS model.
6!! Used to test the inversion routine for a very dense model.
7!! The model is similar to Lincns.f90, but many nonzeros are so small
8!! that they are removed by the factorization routine. The
9!! factorization will therefore have to handle many fairly large
10!! blocks with many pivots and streamlining the updates is a
11!! challenge.
12!!
13!!
14!! For more information about the individual callbacks, please have a look at the source code.
15
16!> Main program. A simple setup and call of CONOPT
17!!
18Program lincns
19 Use proginfop
20 Use coidef
21 Use omp_lib
22 IMPLICIT NONE
23!
24! Declare the user callback routines as Integer, External:
25!
26 Integer, External :: lincns_readmatrix ! Mandatory Matrix definition routine defined below
27 Integer, External :: std_status ! Standard callback for displaying solution status
28 Integer, External :: lincns_solution ! Standard callback for displaying solution values
29 Integer, External :: std_message ! Standard callback for managing messages
30 Integer, External :: std_errmsg ! Standard callback for managing error messages
31 Integer, External :: std_triord ! Standard callback for triangular order
32#if defined(itl)
33!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lincns_ReadMatrix
34!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Status
35!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lincns_Solution
36!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Message
37!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_ErrMsg
38!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_TriOrd
39#endif
40!
41! Control vector
42!
43 INTEGER, Dimension(:), Pointer :: cntvect
44 INTEGER :: coi_error
45!
46! Locals
47!
48 Integer :: n
49
50 real*8 time0, time1, time2
51!
52! Create and initialize a Control Vector
53!
54 call startup
55 coi_error = coi_createfort( cntvect )
56!
57! Define the number of variables and constraints
58!
59 n = 1000
60!
61! Tell CONOPT about the size of the model by populating the Control Vector:
62!
63! Number of variables = N
64!
65 coi_error = max( coi_error, coidef_numvar( cntvect, n ) )
66!
67! Number of equations = M
68!
69 coi_error = max( coi_error, coidef_numcon( cntvect, n ) )
70!
71! Number of nonzeros in the Jacobian: N*N
72!
73 coi_error = max( coi_error, coidef_numnz( cntvect, n*n ) )
74!
75! Number of nonlinear nonzeros. 0 -- the model is linear
76!
77 coi_error = max( coi_error, coidef_numnlnz( cntvect, 0 ) )
78!
79! Square system ( No objective is defined)
80!
81 coi_error = max( coi_error, coidef_square( cntvect, 1 ) )
82!
83! Options file is test.opt
84!
85 coi_error = max( coi_error, coidef_optfile( cntvect, 'mp_lincns2.opt' ) )
86!
87! Tell CONOPT about the callback routines:
88!
89 coi_error = max( coi_error, coidef_readmatrix( cntvect, lincns_readmatrix ) )
90 coi_error = max( coi_error, coidef_status( cntvect, std_status ) )
91 coi_error = max( coi_error, coidef_solution( cntvect, lincns_solution ) )
92 coi_error = max( coi_error, coidef_message( cntvect, std_message ) )
93 coi_error = max( coi_error, coidef_errmsg( cntvect, std_errmsg ) )
94 coi_error = max( coi_error, coidef_triord( cntvect, std_triord ) )
95
96#if defined(LICENSE_INT_1) && defined(LICENSE_INT_2) && defined(LICENSE_INT_3) && defined(LICENSE_TEXT)
97 coi_error = max( coi_error, coidef_license( cntvect, license_int_1, license_int_2, license_int_3, license_text) )
98#endif
99
100 If ( coi_error .ne. 0 ) THEN
101 write(*,*)
102 write(*,*) '**** Fatal Error while loading CONOPT Callback routines.'
103 write(*,*)
104 call flog( "Skipping Solve due to setup errors", 1 )
105 ENDIF
106!
107! Start CONOPT with a single thread:
108!
109 time0 = omp_get_wtime()
110 coi_error = coi_solve( cntvect )
111 time1 = omp_get_wtime() - time0
112
113 write(*,*)
114 write(*,*) 'Single Thread: End of LinCns2. Return code=',coi_error
115
116 If ( coi_error /= 0 ) then
117 call flog( "Single Thread: Errors encountered during solution", 1 )
118 elseif ( stacalls == 0 .or. solcalls == 0 ) then
119 call flog( "Single Thread: Status or Solution routine was not called", 1 )
120 elseif ( mstat /= 15 .or. sstat /= 1 ) then ! 15 since the model is linear
121 call flog( "Single Thread: The model or solver status was not (15,1) as expected", 1 )
122 elseif ( miter /= 1 ) then ! we expect one iteration for a linear CNS
123 call flog( "Single Thread: The iteration count was not 1 as expected", 1 )
124 endif
125!
126! Start CONOPT again using multiple threads everywhere:
127!
128 coi_error = max( coi_error, coidef_threads( cntvect, 0 ) ) ! 0 means use as many threads as you can
129 time0 = omp_get_wtime()
130 coi_error = coi_solve( cntvect )
131 time2 = omp_get_wtime() - time0
132
133 write(*,*)
134 write(*,*) 'Multiple Threads: End of LinCns2. Return code=',coi_error
135
136 If ( coi_error /= 0 ) then
137 call flog( "Multiple Threads: Errors encountered during solution", 1 )
138 elseif ( stacalls == 0 .or. solcalls == 0 ) then
139 call flog( "Multiple Threads: Status or Solution routine was not called", 1 )
140 elseif ( mstat /= 15 .or. sstat /= 1 ) then ! 15 since the model is linear
141 call flog( "Multiple Threads: The model or solver status was not (15,1) as expected", 1 )
142 elseif ( miter /= 1 ) then ! we expect one iteration for a linear CNS
143 call flog( "Multiple Threads: The iteration count was not 1 as expected", 1 )
144 endif
145
146 if ( coi_free( cntvect ) /= 0 ) call flog( "Error while freeing control vector", 1 )
147
148 write(*,*)
149 write(*,"('Time for single thread',f10.3)") time1
150 write(*,"('Time for multi thread',f10.3)") time2
151 write(*,"('Speedup ',f10.3)") time1/time2
152 write(*,"('Efficiency ',f10.3)") time1/time2/omp_get_max_threads()
153
154 call flog( "Successful Solve", 0 )
155
156end Program lincns
157!
158! =====================================================================
159! Define information about the model structure
160!
161
162!> Define information about the model
163!!
164!! @include{doc} readMatrix_params.dox
165Integer Function lincns_readmatrix( lower, curr, upper, vsta, type, rhs, esta, &
166 colsta, rowno, value, nlflag, n, m, nz, usrmem )
167#if defined(itl)
168!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lincns_ReadMatrix
169#endif
170 IMPLICIT NONE
171 integer, intent (in) :: n ! number of variables
172 integer, intent (in) :: m ! number of constraints
173 integer, intent (in) :: nz ! number of nonzeros
174 real*8, intent (in out), dimension(n) :: lower ! vector of lower bounds
175 real*8, intent (in out), dimension(n) :: curr ! vector of initial values
176 real*8, intent (in out), dimension(n) :: upper ! vector of upper bounds
177 integer, intent (in out), dimension(n) :: vsta ! vector of initial variable status
178 ! (not defined here)
179 integer, intent (out), dimension(m) :: type ! vector of equation types
180 integer, intent (in out), dimension(m) :: esta ! vector of initial equation status
181 ! (not defined here)
182 real*8, intent (in out), dimension(m) :: rhs ! vector of right hand sides
183 integer, intent (in out), dimension(n+1) :: colsta ! vector with start of column indices
184 integer, intent (out), dimension(nz) :: rowno ! vector of row numbers
185 integer, intent (in out), dimension(nz) :: nlflag ! vector of nonlinearity flags
186 real*8, intent (in out), dimension(nz) :: value ! vector of matrix values
187 real*8 usrmem(*) ! optional user memory
188
189 Integer :: i, j, k
190!
191! Define the information for the columns. All are unbounded with initial value I
192!
193 do i = 1, n
194 curr(i) = i
195 enddo
196!
197! All are equalities:
198!
199 do i = 1, m
200 type(i) = 0
201 enddo
202!
203! Right hand sides: Initialize to 0 and revise so the solution becomes
204! x(i) = 1.0
205!
206 do i = 1, m
207 rhs(i) = 0.d0
208 enddo
209!
210! Define the structure and content of the Jacobian:
211! We have a pattern with the largest values on the cross-diagonal from
212! (1,N) to (N,1) with gradually smaller value away from this diagonal.
213!
214 k = 1
215 do i = 1, n
216 colsta(i) = k
217 do j = 1, n
218 rowno(k) = j
219 value(k) = 2.0d0**(-abs(i+j-(n+1))) ! many value are really small
220 rhs(j) = rhs(j) + value(k)
221 k = k + 1
222 enddo
223 enddo
224 colsta(n+1) = k
225
227
228end Function lincns_readmatrix
229
230Integer Function lincns_solution( XVAL, XMAR, XBAS, XSTA, YVAL, YMAR, YBAS, YSTA, N, M, USRMEM )
231#if defined(itl)
232!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lincns_Solution
233#endif
234!
235! We write the solution values to the 'Documentation file' on unit 10
236! and test that the solution values are 1.0
237!
238 Use proginfop
239 IMPLICIT NONE
240 INTEGER, Intent(IN) :: n, m
241 INTEGER, Intent(IN), Dimension(N) :: xbas, xsta
242 INTEGER, Intent(IN), Dimension(M) :: ybas, ysta
243 real*8, Intent(IN), Dimension(N) :: xval, xmar
244 real*8, Intent(IN), Dimension(M) :: yval, ymar
245 real*8, Intent(IN OUT) :: usrmem(*)
246
247 INTEGER i
248 CHARACTER*5, Parameter, Dimension(4) :: stat = (/ 'Lower','Upper','Basic','Super' /)
249
250 Logical error
251
252 error = .false.
254 WRITE(10,"(/' Variable Solution value Reduced cost B-stat')")
255 DO i = 1, n
256 WRITE(10,"(1P,I7,E20.6,E16.6,4X,A5 )") i, xval(i), xmar(i), stat(1+xbas(i))
257 If ( abs(xval(i)-1.d0) > 1.d-7 ) then
258 write(10,*) '**** Bad solution value'
259 error = .true.
261 Endif
262 ENDDO
263
264 WRITE(10,"(/' Constrnt Activity level Marginal cost B-stat')")
265 DO i = 1, m
266 WRITE(10,"(1P,I7,E20.6,E16.6,4X,A5 )") i, yval(i), ymar(i), stat(1+ybas(i))
267 ENDDO
268
269 solcalls = solcalls + 1
270
271END Function lincns_solution
272
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_triord(mode, type, status, irow, icol, inf, value, resid, usrmem)
Definition comdecl.f90:291
integer function std_errmsg(rowno, colno, posno, msglen, usrmem, msg)
Definition comdecl.f90:248
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_triord(cntvect, coi_triord)
define callback routine for providing the triangular order information.
integer function coidef_threads(cntvect, threads)
number of threads allowed internally in CONOPT.
integer function coidef_square(cntvect, square)
square models.
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_numnz(cntvect, numnz)
defines the number of nonzero elements in the Jacobian.
Definition coistart.f90:437
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 lincns_readmatrix(lower, curr, upper, vsta, type, rhs, esta, colsta, rowno, value, nlflag, n, m, nz, usrmem)
Define information about the model.
integer function lincns_solution(xval, xmar, xbas, xsta, yval, ymar, ybas, ysta, n, m, usrmem)
program lincns
Main program. A simple setup and call of CONOPT.
Definition mp_lincns.f90:13
subroutine flog(msg, code)
Definition comdeclp.f90:42
integer miter
Definition comdeclp.f90:16
subroutine startup
Definition comdeclp.f90:25