CONOPT
Loading...
Searching...
No Matches
mp_leastsq14.f90
Go to the documentation of this file.
1!> @file mp_leastsq14.f90
2!! @ingroup FORTOPENMP_EXAMPLES
3!!
4!!
5!! This model is an extension of leastsq.f90 in which we test that the
6!! solution with 1, 2, and 4 threads give the same solution when
7!! threadc is defined as 4.
8!!
9!! We solve the following nonlinear least squares model:
10!!
11!! min sum(i, res(i)**2 )
12!!
13!! sum(j, a(i,j)*x(j) + b(i,j)*x(j)**2 ) + res(i) = obs(i)
14!!
15!! where a, b, and obs are known data, and res and x are the
16!! variables of the model.
17!!
18!!
19!! For more information about the individual callbacks, please have a look at the source code.
20
21#if defined(_WIN32) && !defined(_WIN64)
22#define dec_directives_win32
23#endif
24
25module lsq14
26 Integer, parameter :: Nobs = 700
27 Integer, parameter :: DimX = 500
28 real*8, dimension(Nobs,DimX) :: a, b
29 real*8, dimension(Nobs) :: obs
30 real*8, dimension(Nobs+DimX) :: xprim, xprim1, xprim2, xprim4
31end module lsq14
33
34REAL FUNCTION rndx( )
35!
36! Defines a pseudo random number between 0 and 1
37!
38 IMPLICIT NONE
39
40 Integer, save :: seed = 12359
41
42 seed = mod(seed*1027+25,1048576)
43 rndx = float(seed)/float(1048576)
44
45END FUNCTION rndx
46
47subroutine defdata
48!
49! Define values for A, B, and Obs
50!
51 Use lsq14
52 IMPLICIT NONE
53
54 Integer :: i, j
55 real*8, Parameter :: xtarg = -1.0
56 real*8, Parameter :: noise = 1.0
57 Real, External :: Rndx
58 real*8 :: o
59
60 do i = 1, nobs
61 o = 0.d0
62 do j = 1, dimx
63 a(i,j) = rndx()
64 b(i,j) = rndx()
65 o = o + a(i,j) * xtarg + b(i,j) * xtarg**2
66 enddo
67 obs(i) = o + noise * rndx()
68 enddo
69
70end subroutine defdata
71!
72! Main program.
73!
74!> Main program. A simple setup and call of CONOPT
75!!
76Program leastsquare
77
79 Use conopt
80 Use lsq14
81 Use omp_lib
82 implicit None
83!
84! Declare the user callback routines as Integer, External:
85!
86 Integer, External :: lsq_readmatrix ! Mandatory Matrix definition routine defined below
87 Integer, External :: lsq_fdeval ! Function and Derivative evaluation routine
88 ! needed a nonlinear model.
89 Integer, External :: std_status ! Standard callback for displaying solution status
90 Integer, External :: lsq_solution ! Special callback for displaying and testing solution values
91 Integer, External :: std_message ! Standard callback for managing messages
92 Integer, External :: std_errmsg ! Standard callback for managing error messages
93#ifdef dec_directives_win32
94!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_ReadMatrix
95!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_FDEval
96!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Status
97!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_Solution
98!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Message
99!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_ErrMsg
100#endif
101!
102! Control vector
103!
104 INTEGER, Dimension(:), Pointer :: cntvect
105 INTEGER :: coi_error
106
107 real*8 time0, time1, time2, time4
108 Integer :: i, maxthread
109 Logical :: same12, same14, same24
110
111 call startup
112!
113! Define data
114!
115 Call defdata
116!
117! Create and initialize a Control Vector
118!
119 coi_error = coi_create( cntvect )
120!
121! Tell CONOPT about the size of the model by populating the Control Vector:
122!
123 coi_error = max( coi_error, coidef_numvar( cntvect, nobs + dimx ) )
124 coi_error = max( coi_error, coidef_numcon( cntvect, nobs + 1 ) )
125 coi_error = max( coi_error, coidef_numnz( cntvect, nobs * (dimx+2) ) )
126 coi_error = max( coi_error, coidef_numnlnz( cntvect, nobs * (dimx+1) ) )
127 coi_error = max( coi_error, coidef_optdir( cntvect, -1 ) ) ! Minimize
128 coi_error = max( coi_error, coidef_objcon( cntvect, nobs + 1 ) ) ! Objective is last constraint
129 coi_error = max( coi_error, coidef_optfile( cntvect, 'mp_leastsq14.opt' ) )
130!
131! Tell CONOPT about the callback routines:
132!
133 coi_error = max( coi_error, coidef_readmatrix( cntvect, lsq_readmatrix ) )
134 coi_error = max( coi_error, coidef_fdeval( cntvect, lsq_fdeval ) )
135 coi_error = max( coi_error, coidef_status( cntvect, std_status ) )
136 coi_error = max( coi_error, coidef_solution( cntvect, lsq_solution ) )
137 coi_error = max( coi_error, coidef_message( cntvect, std_message ) )
138 coi_error = max( coi_error, coidef_errmsg( cntvect, std_errmsg ) )
139 coi_error = max( coi_error, coidef_debugfv( cntvect, 0 ) ) ! Debug Fdeval on or off
140
141#if defined(CONOPT_LICENSE_INT_1) && defined(CONOPT_LICENSE_INT_2) && defined(CONOPT_LICENSE_INT_3) && defined(CONOPT_LICENSE_TEXT)
142 coi_error = max( coi_error, coidef_license( cntvect, conopt_license_int_1, conopt_license_int_2, conopt_license_int_3, conopt_license_text) )
143#endif
144
145 If ( coi_error .ne. 0 ) THEN
146 write(*,*)
147 write(*,*) '**** Fatal Error while loading CONOPT Callback routines.'
148 write(*,*)
149 call flog( "Skipping Solve due to setup errors", 1 )
150 ENDIF
151 maxthread = omp_get_max_threads() ! Define ThreadC as the largest number of threads we can get
152 coi_error = max( coi_error, coidef_threadc( cntvect, maxthread ) )
153!
154! Start CONOPT -- First time using a single thread:
155!
156 coi_error = max( coi_error, coidef_threads( cntvect, 1 ) ) ! 1 means use 1 thread
157 time0 = omp_get_wtime()
158 coi_error = coi_solve( cntvect )
159 time1 = omp_get_wtime() - time0
160
161 write(*,*)
162 write(*,*) 'End of Least Square example 1 with 1 thread. Return code=',coi_error
163
164 If ( coi_error /= 0 ) then
165 call flog( "One thread: Errors encountered during solution", 1 )
166 elseif ( stacalls == 0 .or. solcalls == 0 ) then
167 call flog( "One thread: Status or Solution routine was not called", 1 )
168 elseif ( .not. ( sstat == 1 .and. mstat == 2 ) ) then
169 call flog( "One thread: Solver or Model status was not as expected (1,2)", 1 )
170 elseif ( abs( obj - 19.44434311d0 ) > 1.d-7 ) then
171 call flog( "One thread: Incorrect objective returned", 1 )
172 endif
173 xprim1 = xprim
174!
175! Start CONOPT again -- this time using 2 threads:
176!
177 coi_error = max( coi_error, coidef_threads( cntvect, 2 ) ) ! 2 means use 2 threads
178 time0 = omp_get_wtime()
179 coi_error = coi_solve( cntvect )
180 time2 = omp_get_wtime() - time0
181
182 write(*,*)
183 write(*,*) 'Two threads: End of Least Square example 1. Return code=',coi_error
184
185 If ( coi_error /= 0 ) then
186 call flog( "Two threads: Errors encountered during solution", 1 )
187 elseif ( stacalls == 0 .or. solcalls == 0 ) then
188 call flog( "Two threads: Status or Solution routine was not called", 1 )
189 elseif ( .not. ( sstat == 1 .and. mstat == 2 ) ) then
190 call flog( "Two threads: Solver or Model status was not as expected (1,2)", 1 )
191 elseif ( abs( obj - 19.44434311d0 ) > 1.d-7 ) then
192 call flog( "Two threads: Incorrect objective returned", 1 )
193 endif
194 xprim2 = xprim
195 if ( maxthread >= 4 ) then
196!
197! Start CONOPT again -- this time using 4 threads (if available on this machine):
198!
199 coi_error = max( coi_error, coidef_threads( cntvect, 4 ) ) ! 4 means use 4 threads
200 time0 = omp_get_wtime()
201 coi_error = coi_solve( cntvect )
202 time4 = omp_get_wtime() - time0
203
204 write(*,*)
205 write(*,*) 'Four threads: End of Least Square example 1. Return code=',coi_error
206
207 If ( coi_error /= 0 ) then
208 call flog( "Four threads: Errors encountered during solution", 1 )
209 elseif ( stacalls == 0 .or. solcalls == 0 ) then
210 call flog( "Four threads: Status or Solution routine was not called", 1 )
211 elseif ( .not. ( sstat == 1 .and. mstat == 2 ) ) then
212 call flog( "Four threads: Solver or Model status was not as expected (1,2)", 1 )
213 elseif ( abs( obj - 19.44434311d0 ) > 1.d-7 ) then
214 call flog( "Four threads: Incorrect objective returned", 1 )
215 endif
216 xprim4 = xprim
217 Endif
218
219 same12 = .true.
220 same14 = .true.
221 same24 = .true.
222 Do i = 1, nobs+dimx
223 If ( xprim2(i) /= xprim1(i) ) then
224 if ( same12 ) then
225 write(*,*) 'Solution 2 is different from Solution 1. First diff for index',i
226 endif
227 same12 = .false.
228 endif
229 if ( maxthread >= 4 ) then
230 If ( xprim4(i) /= xprim1(i) ) then
231 if ( same14 ) then
232 write(*,*) 'Solution 4 is different from Solution 1. First diff for index',i
233 endif
234 same14 = .false.
235 endif
236 If ( xprim4(i) /= xprim2(i) ) then
237 if ( same24 ) then
238 write(*,*) 'Solution 4 is different from Solution 2. First diff for index',i
239 endif
240 same24 = .false.
241 endif
242 endif
243 Enddo
244 write(*,*)
245 write(*,"('Time for single thread ',f10.3)") time1
246 write(*,"('Time for two threads',f10.3)") time2
247 if ( maxthread >= 4 ) &
248 write(*,"('Time for four threads',f10.3)") time4
249 write(*,"('Speedup two threads',f10.3)") time1/time2
250 write(*,"('Efficiency two threads',f10.3)") time1/time2/2
251 if ( maxthread >= 4 ) &
252 write(*,"('Speedup four threads',f10.3)") time1/time4
253 if ( maxthread >= 4 ) &
254 write(*,"('Efficiency four threads',f10.3)") time1/time4/4
255
256#if defined (nocomp)
257#else
258 if ( same12 .and. same14 .and. same24 ) then
259 write(*,*) 'All solutions are the same'
260 else
261 write(*,*) 'The solutions are NOT the same'
262 call flog( "Solutions are not the same", 1 )
263 Endif
264#endif
265
266 if ( coi_free( cntvect ) /= 0 ) call flog( "Error while freeing control vector", 1 )
268 call flog( "Successful Solve", 0 )
269
270End Program leastsquare
271!
272! ============================================================================
273! Define information about the model:
274!
275
276!> Define information about the model
277!!
278!! @include{doc} readMatrix_params.dox
279Integer Function lsq_readmatrix( lower, curr, upper, vsta, type, rhs, esta, &
280 colsta, rowno, value, nlflag, n, m, nz, &
281 usrmem )
282#ifdef dec_directives_win32
283!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_ReadMatrix
284#endif
285 Use lsq14
286 implicit none
287 integer, intent (in) :: n ! number of variables
288 integer, intent (in) :: m ! number of constraints
289 integer, intent (in) :: nz ! number of nonzeros
290 real*8, intent (in out), dimension(n) :: lower ! vector of lower bounds
291 real*8, intent (in out), dimension(n) :: curr ! vector of initial values
292 real*8, intent (in out), dimension(n) :: upper ! vector of upper bounds
293 integer, intent (in out), dimension(n) :: vsta ! vector of initial variable status
294 ! (not defined here)
295 integer, intent (out), dimension(m) :: type ! vector of equation types
296 integer, intent (in out), dimension(m) :: esta ! vector of initial equation status
297 ! (not defined here)
298 real*8, intent (in out), dimension(m) :: rhs ! vector of right hand sides
299 integer, intent (in out), dimension(n+1) :: colsta ! vector with start of column indices
300 integer, intent (out), dimension(nz) :: rowno ! vector of row numbers
301 integer, intent (in out), dimension(nz) :: nlflag ! vector of nonlinearity flags
302 real*8, intent (in out), dimension(nz) :: value ! vector of matrix values
303 real*8 usrmem(*) ! optional user memory
304
305 Integer :: i, j, k
306!
307! Information about Variables:
308! Default: Lower = -Inf, Curr = 0, and Upper = +inf.
309! Default: the status information in Vsta is not used.
310!
311 do i = 1, dimx
312 curr(i) = -0.8d0
313 enddo
314!
315! Information about Constraints:
316! Default: Rhs = 0
317! Default: the status information in Esta and the function
318! value in FV are not used.
319! Default: Type: There is no default.
320! 0 = Equality,
321! 1 = Greater than or equal,
322! 2 = Less than or equal,
323! 3 = Non binding.
324!
325! Constraints 1 to Nobs:
326! Rhs = Obs(i) and type Equality
327!
328 do i = 1, nobs
329 rhs(i) = obs(i)
330 type(i) = 0
331 enddo
332!
333! Constraint Nobs + 1 (Objective)
334! Rhs = 0 and type Non binding
335!
336 type(nobs+1) = 3
337!
338! Information about the Jacobian. CONOPT expects a columnwise
339! representation in Rowno, Value, Nlflag and Colsta.
340!
341! Colsta = Start of column indices (No Defaults):
342! Rowno = Row indices
343! Value = Value of derivative (by default only linear
344! derivatives are used)
345! Nlflag = 0 for linear and 1 for nonlinear derivative
346! (not needed for completely linear models)
347!
348!
349! Indices
350! x(i) res(j)
351! j: NL L=1
352! obj: NL
353!
354 k = 1
355 do j = 1, dimx
356 colsta(j) = k
357 do i = 1, nobs
358 rowno(k) = i
359 nlflag(k) = 1
360 k = k + 1
361 enddo
362 enddo
363 do i = 1, nobs
364 colsta(dimx+i) = k
365 rowno(k) = i
366 nlflag(k) = 0
367 value(k) = 1.d0
368 k = k + 1
369 rowno(k) = nobs+1
370 nlflag(k) = 1
371 k = k + 1
372 enddo
373 colsta(dimx+nobs+1) = k
374
375 lsq_readmatrix = 0 ! Return value means OK
376
377end Function lsq_readmatrix
378!
379!==========================================================================
380! Compute nonlinear terms and non-constant Jacobian elements
381!
382
383!> Compute nonlinear terms and non-constant Jacobian elements
384!!
385!! @include{doc} fdeval_params.dox
386Integer Function lsq_fdeval( x, g, jac, rowno, jcnm, mode, ignerr, errcnt, &
387 n, nz, thread, usrmem )
388#ifdef dec_directives_win32
389!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_FDEval
390#endif
391 use lsq14
392 implicit none
393 integer, intent (in) :: n ! number of variables
394 integer, intent (in) :: rowno ! number of the row to be evaluated
395 integer, intent (in) :: nz ! number of nonzeros in this row
396 real*8, intent (in), dimension(n) :: x ! vector of current solution values
397 real*8, intent (in out) :: g ! constraint value
398 real*8, intent (in out), dimension(n) :: jac ! vector of derivatives for current constraint
399 integer, intent (in), dimension(nz) :: jcnm ! list of variables that appear nonlinearly
400 ! in this row. Ffor information only.
401 integer, intent (in) :: mode ! evaluation mode: 1 = function value
402 ! 2 = derivatives, 3 = both
403 integer, intent (in) :: ignerr ! if 1 then errors can be ignored as long
404 ! as errcnt is incremented
405 integer, intent (in out) :: errcnt ! error counter to be incremented in case
406 ! of function evaluation errors.
407 integer, intent (in) :: thread
408 real*8 usrmem(*) ! optional user memory
409
410 integer :: i,j
411 real*8 :: s
412!
413! Row Nobs+1: the objective function is nonlinear
414!
415 if ( rowno .eq. nobs+1 ) then
416!
417! Mode = 1 or 3. Function value: G = sum(i, res(i)**2)
418!
419 if ( mode .eq. 1 .or. mode .eq. 3 ) then
420 s = 0.d0
421 do i = 1, nobs
422 s = s + x(dimx+i)**2
423 enddo
424 g = s
425 endif
426!
427! Mode = 2 or 3: Derivative values:
428!
429 if ( mode .eq. 2 .or. mode .eq. 3 ) then
430 do i = 1, nobs
431 jac(dimx+i) = 2*x(dimx+i)
432 enddo
433 endif
434!
435! Row 1 to Nobs: The observation definitions:
436!
437 else
438!
439! Mode = 1 or 3: Function value - x-part only
440!
441 if ( mode .eq. 1 .or. mode .eq. 3 ) then
442 s = 0.d0
443 do j = 1, dimx
444 s = s + a(rowno,j)*x(j) + b(rowno,j)*x(j)**2
445 enddo
446 g = s
447 endif
448!
449! Mode = 2 or 3: Derivatives
450!
451 if ( mode .eq. 2 .or. mode .eq. 3 ) then
452 do j = 1, dimx
453 jac(j) = a(rowno,j) + 2.d0*b(rowno,j)*x(j)
454 enddo
455 endif
456 endif
457 lsq_fdeval = 0
458
459end Function lsq_fdeval
460
461Integer Function lsq_solution( XVAL, XMAR, XBAS, XSTA, YVAL, YMAR, YBAS, YSTA, N, M, USRMEM )
462#ifdef dec_directives_win32
463!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_Solution
464#endif
465!
466! Simple implementation in which we write the solution values to
467! the 'Documentation file' on unit 10.
468!
469 Use proginfop
470 Use lsq14
471 IMPLICIT NONE
472 INTEGER, Intent(IN) :: n, m
473 INTEGER, Intent(IN), Dimension(N) :: xbas, xsta
474 INTEGER, Intent(IN), Dimension(M) :: ybas, ysta
475 real*8, Intent(IN), Dimension(N) :: xval, xmar
476 real*8, Intent(IN), Dimension(M) :: yval, ymar
477 real*8, Intent(IN OUT) :: usrmem(*)
478
479 INTEGER i
480 CHARACTER*5, Parameter, Dimension(4) :: stat = (/ 'Lower','Upper','Basic','Super' /)
481
482 WRITE(10,"(/' Variable Solution value Reduced cost B-stat')")
483 DO i = 1, n
484 WRITE(10,"(1P,I7,E20.6,E16.6,4X,A5 )") i, xval(i), xmar(i), stat(1+xbas(i))
485 ENDDO
486
487 WRITE(10,"(/' Constrnt Activity level Marginal cost B-stat')")
488 DO i = 1, m
489 WRITE(10,"(1P,I7,E20.6,E16.6,4X,A5 )") i, yval(i), ymar(i), stat(1+ybas(i))
490 ENDDO
491 xprim(1:n) = xval(1:n)
492 Call flush(10)
493
494 solcalls = solcalls + 1
495 lsq_solution = 0
496
497END Function lsq_solution
integer function std_status(modsta, solsta, iter, objval, usrmem)
Definition comdecl.f90:88
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_threadc(cntvect, threadc)
check for thread compatibility.
Definition conopt.f90:727
integer(c_int) function coidef_debugfv(cntvect, debugfv)
turn Debugging of FDEval on and off.
Definition conopt.f90:387
integer(c_int) function coidef_threads(cntvect, threads)
number of threads allowed internally in CONOPT.
Definition conopt.f90:627
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
integer function lsq_readmatrix(lower, curr, upper, vsta, type, rhs, esta, colsta, rowno, value, nlflag, n, m, nz, usrmem)
Define information about the model.
Definition leastl1.f90:166
integer function lsq_fdeval(x, g, jac, rowno, jcnm, mode, ignerr, errcnt, n, nz, thread, usrmem)
Compute nonlinear terms and non-constant Jacobian elements.
Definition leastl1.f90:287
#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 leastsquare
Main program. A simple setup and call of CONOPT.
Definition leastsq.f90:68
integer function lsq_solution(xval, xmar, xbas, xsta, yval, ymar, ybas, ysta, n, m, usrmem)
real *8, dimension(nobs+dimx) xprim1
real *8, dimension(nobs+dimx) xprim4
real *8, dimension(nobs, dimx) b
real *8, dimension(nobs+dimx) xprim
real *8, dimension(nobs+dimx) xprim2
real *8, dimension(nobs, dimx) a
real *8, dimension(nobs) obs
subroutine flog(msg, code)
Definition comdeclp.f90:48
real *8 obj
Definition comdeclp.f90:19
subroutine startup
Definition comdeclp.f90:31