CONOPT
Loading...
Searching...
No Matches
mp_leastsq5.f90
Go to the documentation of this file.
1!> @file mp_leastsq5.f90
2!! @ingroup FORTOPENMP_EXAMPLES
3!!
4!! @copydoc leastsq5.f90
5
6#if defined(_WIN32) && !defined(_WIN64)
7#define dec_directives_win32
8#endif
9
10module lsq5
11 Integer, parameter :: Nobs = 700
12 Integer, parameter :: DimX = 500
13 real*8, dimension(Nobs,DimX) :: a, b
14 real*8, dimension(Nobs) :: obs
15end module lsq5
17
18REAL FUNCTION rndx( )
19!
20! Defines a pseudo random number between 0 and 1
21!
22 IMPLICIT NONE
23
24 Integer, save :: seed = 12359
25
26 seed = mod(seed*1027+25,1048576)
27 rndx = float(seed)/float(1048576)
28
29END FUNCTION rndx
30
31subroutine defdata
32!
33! Define values for A, B, and Obs
34!
35 Use lsq5
36 IMPLICIT NONE
37
38 Integer :: i, j
39 real*8, Parameter :: xtarg = -1.0
40 real*8, Parameter :: noise = 1.0
41 Real, External :: Rndx
42 real*8 :: o
43
44 do i = 1, nobs
45 o = 0.d0
46 do j = 1, dimx
47 a(i,j) = rndx()
48 b(i,j) = rndx()
49 o = o + a(i,j) * xtarg + b(i,j) * xtarg**2
50 enddo
51 obs(i) = o + noise * rndx()
52 enddo
53
54end subroutine defdata
55!
56! Main program.
57!
58!> Main program. A simple setup and call of CONOPT
59!!
60Program leastsquare
61
63 Use conopt
64 Use lsq5
65 Use omp_lib
66 implicit None
67!
68! Declare the user callback routines as Integer, External:
69!
70 Integer, External :: lsq_readmatrix ! Mandatory Matrix definition routine defined below
71 Integer, External :: lsq_fdeval ! Function and Derivative evaluation routine
72 ! needed a nonlinear model.
73 Integer, External :: lsq_2ddir ! Directional 2nd derivative routine
74 Integer, External :: std_status ! Standard callback for displaying solution status
75 Integer, External :: std_solution ! Standard callback for displaying solution values
76 Integer, External :: std_message ! Standard callback for managing messages
77 Integer, External :: std_errmsg ! Standard callback for managing error messages
78#ifdef dec_directives_win32
79!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_ReadMatrix
80!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_FDEval
81!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_2DDir
82!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Status
83!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Solution
84!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Message
85!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_ErrMsg
86#endif
87!
88! Control vector
89!
90 INTEGER, Dimension(:), Pointer :: cntvect
91 INTEGER :: coi_error
92
93 real*8 time0, time1, time2, time3
94
95 call startup
96!
97! Define data
98!
99 Call defdata
100!
101! Create and initialize a Control Vector
102!
103 coi_error = coi_create( cntvect )
104!
105! Tell CONOPT about the size of the model by populating the Control Vector:
106!
107 coi_error = max( coi_error, coidef_numvar( cntvect, nobs + dimx ) )
108 coi_error = max( coi_error, coidef_numcon( cntvect, nobs + 1 ) )
109 coi_error = max( coi_error, coidef_numnz( cntvect, nobs * (dimx+2) ) )
110 coi_error = max( coi_error, coidef_numnlnz( cntvect, nobs * (dimx+1) ) )
111 coi_error = max( coi_error, coidef_optdir( cntvect, -1 ) ) ! Minimize
112 coi_error = max( coi_error, coidef_objcon( cntvect, nobs + 1 ) ) ! Objective is last constraint
113 coi_error = max( coi_error, coidef_optfile( cntvect, 'mp_leastsq5.opt' ) )
114!
115! Tell CONOPT about the callback routines:
116!
117 coi_error = max( coi_error, coidef_readmatrix( cntvect, lsq_readmatrix ) )
118 coi_error = max( coi_error, coidef_fdeval( cntvect, lsq_fdeval ) )
119 coi_error = max( coi_error, coidef_2ddir( cntvect, lsq_2ddir ) )
120 coi_error = max( coi_error, coidef_status( cntvect, std_status ) )
121 coi_error = max( coi_error, coidef_solution( cntvect, std_solution ) )
122 coi_error = max( coi_error, coidef_message( cntvect, std_message ) )
123 coi_error = max( coi_error, coidef_errmsg( cntvect, std_errmsg ) )
124 coi_error = max( coi_error, coidef_debugfv( cntvect, 0 ) ) ! Debug Fdeval on or off
125
126#if defined(CONOPT_LICENSE_INT_1) && defined(CONOPT_LICENSE_INT_2) && defined(CONOPT_LICENSE_INT_3) && defined(CONOPT_LICENSE_TEXT)
127 coi_error = max( coi_error, coidef_license( cntvect, conopt_license_int_1, conopt_license_int_2, conopt_license_int_3, conopt_license_text) )
128#endif
129
130 If ( coi_error .ne. 0 ) THEN
131 write(*,*)
132 write(*,*) '**** Fatal Error while loading CONOPT Callback routines.'
133 write(*,*)
134 call flog( "Skipping Solve due to setup errors", 1 )
135 ENDIF
136!
137! Start CONOPT with a single thread:
138!
139 time0 = omp_get_wtime()
140 coi_error = coi_solve( cntvect )
141 time1 = omp_get_wtime() - time0
142
143 write(*,*)
144 write(*,*) 'Single Thread: End of Least Square example 5. Return code=',coi_error
145
146 If ( coi_error /= 0 ) then
147 call flog( "Single Thread: Errors encountered during solution", 1 )
148 elseif ( stacalls == 0 .or. solcalls == 0 ) then
149 call flog( "Single Thread: Status or Solution routine was not called", 1 )
150 elseif ( sstat /= 1 .or. mstat /= 2 ) then
151 call flog( "Single Thread: Solver and Model Status was not as expected (1,2)", 1 )
152 elseif ( abs( obj - 19.44434311d0 ) > 1.d-7 ) then
153 call flog( "Single Thread: Incorrect objective returned", 1 )
154 endif
155!
156! Start CONOPT again this time with multiple threads:
157!
158 coi_error = max( coi_error, coidef_threads( cntvect, 0 ) ) ! 0 means use as many threads as you can
159 time0 = omp_get_wtime()
160 coi_error = coi_solve( cntvect )
161 time2 = omp_get_wtime() - time0
162
163 write(*,*)
164 write(*,*) 'Multiple Threads: End of Least Square example 5. Return code=',coi_error
165
166 If ( coi_error /= 0 ) then
167 call flog( "Multiple Threads: Errors encountered during solution", 1 )
168 elseif ( stacalls == 0 .or. solcalls == 0 ) then
169 call flog( "Multiple Threads: Status or Solution routine was not called", 1 )
170 elseif ( sstat /= 1 .or. mstat /= 2 ) then
171 call flog( "Multiple Threads: Solver and Model Status was not as expected (1,2)", 1 )
172 elseif ( abs( obj - 19.44434311d0 ) > 1.d-7 ) then
173 call flog( "Multiple Threads: Incorrect objective returned", 1 )
174 endif
175!
176! Start CONOPT again this time with multiple threads except for 2DDir:
177!
178 coi_error = max( coi_error, coidef_thread2d( cntvect, 1 ) ) ! 1 means use only 1 thread for 2DDir
179 time0 = omp_get_wtime()
180 coi_error = coi_solve( cntvect )
181 time3 = omp_get_wtime() - time0
182
183 write(*,*)
184 write(*,*) 'Multiple Threads -2DDir: End of Least Square example 5. Return code=',coi_error
185
186 If ( coi_error /= 0 ) then
187 call flog( "Multiple Threads -2DDir: Errors encountered during solution", 1 )
188 elseif ( stacalls == 0 .or. solcalls == 0 ) then
189 call flog( "Multiple Threads -2DDir: Status or Solution routine was not called", 1 )
190 elseif ( sstat /= 1 .or. mstat /= 2 ) then
191 call flog( "Multiple Threads -2DDir: Solver and Model Status was not as expected (1,2)", 1 )
192 elseif ( abs( obj - 19.44434311d0 ) > 1.d-7 ) then
193 call flog( "Multiple Threads -2DDir: Incorrect objective returned", 1 )
194 endif
195
196 if ( coi_free( cntvect ) /= 0 ) call flog( "Error while freeing control vector", 1 )
197
198 write(*,*)
199 write(*,"('Time for single thread ',f10.3)") time1
200 write(*,"('Time for multi thread incl 2DDir',f10.3)") time2
201 write(*,"('Time for multi thread excl 2DDir',f10.3)") time3
202 write(*,"('Speedup incl 2DDir',f10.3)") time1/time2
203 write(*,"('Efficiency incl 2DDir',f10.3)") time1/time2/omp_get_max_threads()
204 write(*,"('Speedup excl 2DDir',f10.3)") time1/time3
205 write(*,"('Efficiency excl 2DDir',f10.3)") time1/time3/omp_get_max_threads()
206
207 call flog( "Successful Solve", 0 )
209End Program leastsquare
210!
211! ============================================================================
212! Define information about the model:
213!
214
215!> Define information about the model
216!!
217!! @include{doc} readMatrix_params.dox
218Integer Function lsq_readmatrix( lower, curr, upper, vsta, type, rhs, esta, &
219 colsta, rowno, value, nlflag, n, m, nz, &
220 usrmem )
221#ifdef dec_directives_win32
222!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_ReadMatrix
223#endif
224 Use lsq5
225 implicit none
226 integer, intent (in) :: n ! number of variables
227 integer, intent (in) :: m ! number of constraints
228 integer, intent (in) :: nz ! number of nonzeros
229 real*8, intent (in out), dimension(n) :: lower ! vector of lower bounds
230 real*8, intent (in out), dimension(n) :: curr ! vector of initial values
231 real*8, intent (in out), dimension(n) :: upper ! vector of upper bounds
232 integer, intent (in out), dimension(n) :: vsta ! vector of initial variable status
233 ! (not defined here)
234 integer, intent (out), dimension(m) :: type ! vector of equation types
235 integer, intent (in out), dimension(m) :: esta ! vector of initial equation status
236 ! (not defined here)
237 real*8, intent (in out), dimension(m) :: rhs ! vector of right hand sides
238 integer, intent (in out), dimension(n+1) :: colsta ! vector with start of column indices
239 integer, intent (out), dimension(nz) :: rowno ! vector of row numbers
240 integer, intent (in out), dimension(nz) :: nlflag ! vector of nonlinearity flags
241 real*8, intent (in out), dimension(nz) :: value ! vector of matrix values
242 real*8 usrmem(*) ! optional user memory
243
244 Integer :: i, j, k
245!
246! Information about Variables:
247! Default: Lower = -Inf, Curr = 0, and Upper = +inf.
248! Default: the status information in Vsta is not used.
249!
250 do i = 1, dimx
251 curr(i) = -0.8d0
252 enddo
253!
254! Information about Constraints:
255! Default: Rhs = 0
256! Default: the status information in Esta and the function
257! value in FV are not used.
258! Default: Type: There is no default.
259! 0 = Equality,
260! 1 = Greater than or equal,
261! 2 = Less than or equal,
262! 3 = Non binding.
263!
264! Constraints 1 to Nobs:
265! Rhs = Obs(i) and type Equality
266!
267 do i = 1, nobs
268 rhs(i) = obs(i)
269 type(i) = 0
270 enddo
271!
272! Constraint Nobs + 1 (Objective)
273! Rhs = 0 and type Non binding
274!
275 type(nobs+1) = 3
276!
277! Information about the Jacobian. CONOPT expects a columnwise
278! representation in Rowno, Value, Nlflag and Colsta.
279!
280! Colsta = Start of column indices (No Defaults):
281! Rowno = Row indices
282! Value = Value of derivative (by default only linear
283! derivatives are used)
284! Nlflag = 0 for linear and 1 for nonlinear derivative
285! (not needed for completely linear models)
286!
287!
288! Indices
289! x(i) res(j)
290! j: NL L=1
291! obj: NL
292! 3:
293!
294 k = 1
295 do j = 1, dimx
296 colsta(j) = k
297 do i = 1, nobs
298 rowno(k) = i
299 nlflag(k) = 1
300 k = k + 1
301 enddo
302 enddo
303 do i = 1, nobs
304 colsta(dimx+i) = k
305 rowno(k) = i
306 nlflag(k) = 0
307 value(k) = 1.d0
308 k = k + 1
309 rowno(k) = nobs+1
310 nlflag(k) = 1
311 k = k + 1
312 enddo
313 colsta(dimx+nobs+1) = k
314
315 lsq_readmatrix = 0 ! Return value means OK
316
317end Function lsq_readmatrix
318!
319!==========================================================================
320! Compute nonlinear terms and non-constant Jacobian elements
321!
322
323!> Compute nonlinear terms and non-constant Jacobian elements
324!!
325!! @include{doc} fdeval_params.dox
326Integer Function lsq_fdeval( x, g, jac, rowno, jcnm, mode, ignerr, errcnt, &
327 n, nz, thread, usrmem )
328#ifdef dec_directives_win32
329!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_FDEval
330#endif
331 use lsq5
332 implicit none
333 integer, intent (in) :: n ! number of variables
334 integer, intent (in) :: rowno ! number of the row to be evaluated
335 integer, intent (in) :: nz ! number of nonzeros in this row
336 real*8, intent (in), dimension(n) :: x ! vector of current solution values
337 real*8, intent (in out) :: g ! constraint value
338 real*8, intent (in out), dimension(n) :: jac ! vector of derivatives for current constraint
339 integer, intent (in), dimension(nz) :: jcnm ! list of variables that appear nonlinearly
340 ! in this row. Ffor information only.
341 integer, intent (in) :: mode ! evaluation mode: 1 = function value
342 ! 2 = derivatives, 3 = both
343 integer, intent (in) :: ignerr ! if 1 then errors can be ignored as long
344 ! as errcnt is incremented
345 integer, intent (in out) :: errcnt ! error counter to be incremented in case
346 ! of function evaluation errors.
347 integer, intent (in) :: thread
348 real*8 usrmem(*) ! optional user memory
349
350 integer :: i,j
351 real*8 :: s
352!
353! Row 1: the objective function is nonlinear
354!
355 if ( rowno .eq. nobs+1 ) then
356!
357! Mode = 1 or 3. Function value: G = sum(i, res(i)**2)
358!
359 if ( mode .eq. 1 .or. mode .eq. 3 ) then
360 s = 0.d0
361 do i = 1, nobs
362 s = s + x(dimx+i)**2
363 enddo
364 g = s
365 endif
366!
367! Mode = 2 or 3: Derivative values:
368!
369 if ( mode .eq. 2 .or. mode .eq. 3 ) then
370 do i = 1, nobs
371 jac(dimx+i) = 2*x(dimx+i)
372 enddo
373 endif
374!
375! Row 1 to Nobs: The observation definitions:
376!
377 else
378!
379! Mode = 1 or 3: Function value - x-part only
380!
381 if ( mode .eq. 1 .or. mode .eq. 3 ) then
382 s = 0.d0
383 do j = 1, dimx
384 s = s + a(rowno,j)*x(j) + b(rowno,j)*x(j)**2
385 enddo
386 g = s
387 endif
388!
389! Mode = 2 or 3: Derivatives
390!
391 if ( mode .eq. 2 .or. mode .eq. 3 ) then
392 do j = 1, dimx
393 jac(j) = a(rowno,j) + 2.d0*b(rowno,j)*x(j)
394 enddo
395 endif
396 endif
397 lsq_fdeval = 0
398
399end Function lsq_fdeval
400
401
402!> Computes the second derivative of a constraint in a direction
403!!
404!! @include{doc} 2DDir_params.dox
405Integer Function lsq_2ddir( X, DX, D2G, ROWNO, JCNM, NODRV, N, NJ, THREAD, USRMEM )
406#ifdef dec_directives_win32
407!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_2DDir
408#endif
409 Use lsq5
410 implicit none
411 INTEGER, Intent(IN) :: rowno, n, nj, thread
412 INTEGER, Intent(IN OUT) :: nodrv
413 INTEGER, Dimension(NJ), Intent(IN) :: jcnm
414 real*8, Dimension(N), Intent(IN) :: x
415 real*8, Dimension(N), Intent(IN) :: dx
416 real*8, Dimension(N), Intent(OUT) :: d2g
417 real*8, Intent(IN OUT) :: usrmem(*)
418 integer :: i, j
419
420 if ( rowno .eq. nobs+1 ) then
421!
422! Objective row: sum(i, res(i)**2 )
423!
424 do j = 1, dimx
425 d2g(j) = 0.d0
426 enddo
427 do i = 1, nobs
428 d2g(dimx+i) = 2.d0*dx(dimx+i)
429 enddo
430 else
431!
432! Constraint row with b(i,j)*x(j)**2 as only nonlinear term
433!
434 do j = 1, dimx
435 d2g(j) = 2.d0*b(rowno,j)*dx(j)
436 enddo
437 do i = 1, nobs
438 d2g(dimx+i) = 0.d0
439 enddo
440 endif
441
442 lsq_2ddir = 0
443
444End Function lsq_2ddir
integer function std_solution(xval, xmar, xbas, xsta, yval, ymar, ybas, ysta, n, m, usrmem)
Definition comdecl.f90:132
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_2ddir(cntvect, coi_2ddir)
define callback routine for computing the second derivative for a constraint in a direction.
Definition conopt.f90:1421
integer(c_int) function coidef_license(cntvect, licint1, licint2, licint3, licstring)
define the License Information.
Definition conopt.f90:293
integer(c_int) function coidef_thread2d(cntvect, thread2d)
number of threads allowed for simultaneous 2DDir calls.
Definition conopt.f90:690
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
integer function lsq_2ddir(x, dx, d2g, rowno, jcnm, nodrv, n, nj, thread, usrmem)
Computes the second derivative of a constraint in a direction.
#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
#define nj
Definition mp_trans.c:46
real *8, dimension(nobs) obs
real *8, dimension(nobs, dimx) a
real *8, dimension(nobs, dimx) b
subroutine flog(msg, code)
Definition comdeclp.f90:48
real *8 obj
Definition comdeclp.f90:19
subroutine startup
Definition comdeclp.f90:31