CONOPT
Loading...
Searching...
No Matches
mp_leastsq2.f90
Go to the documentation of this file.
1!> @file mp_leastsq2.f90
2!! @ingroup FORTOPENMP_EXAMPLES
3!!
4!! @copydoc leastsq2.f90
5
6#if defined(_WIN32) && !defined(_WIN64)
7#define dec_directives_win32
8#endif
9
10module lsq2
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 lsq2
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 lsq2
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 lsq2
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 :: std_status ! Standard callback for displaying solution status
74 Integer, External :: std_solution ! Standard callback for displaying solution values
75 Integer, External :: std_message ! Standard callback for managing messages
76 Integer, External :: std_errmsg ! Standard callback for managing error messages
77 Integer, External :: lsq_2dlagrstr ! Callback for second derivatives structure
78 Integer, External :: lsq_2dlagrval ! Callback for second derivatives values
79#ifdef dec_directives_win32
80!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_ReadMatrix
81!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_FDEval
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!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_2DLagrStr
87!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_2DLagrVal
88#endif
89!
90! Control vector
91!
92 INTEGER, Dimension(:), Pointer :: cntvect
93 INTEGER :: coi_error
94
95 real*8 time0, time1, time2
96
97 call startup
98!
99! Define data
100!
101 Call defdata
102!
103! Create and initialize a Control Vector
104!
105 coi_error = coi_create( cntvect )
106!
107! Tell CONOPT about the size of the model by populating the Control Vector:
108!
109 coi_error = max( coi_error, coidef_numvar( cntvect, nobs + dimx ) )
110 coi_error = max( coi_error, coidef_numcon( cntvect, nobs + 1 ) )
111 coi_error = max( coi_error, coidef_numnz( cntvect, nobs * (dimx+2) ) )
112 coi_error = max( coi_error, coidef_numnlnz( cntvect, nobs * (dimx+1) ) )
113 coi_error = max( coi_error, coidef_numhess( cntvect, nobs + dimx ) )
114 coi_error = max( coi_error, coidef_optdir( cntvect, -1 ) ) ! Minimize
115 coi_error = max( coi_error, coidef_objcon( cntvect, nobs + 1 ) ) ! Objective is last constraint
116 coi_error = max( coi_error, coidef_debugfv( cntvect, 0 ) ) ! Debug Fdeval on or off
117 coi_error = max( coi_error, coidef_optfile( cntvect, 'mp_lestsq2.opt' ) )
118!
119! Tell CONOPT about the callback routines:
120!
121 coi_error = max( coi_error, coidef_readmatrix( cntvect, lsq_readmatrix ) )
122 coi_error = max( coi_error, coidef_fdeval( cntvect, lsq_fdeval ) )
123 coi_error = max( coi_error, coidef_status( cntvect, std_status ) )
124 coi_error = max( coi_error, coidef_solution( cntvect, std_solution ) )
125 coi_error = max( coi_error, coidef_message( cntvect, std_message ) )
126 coi_error = max( coi_error, coidef_errmsg( cntvect, std_errmsg ) )
127 coi_error = max( coi_error, coidef_2dlagrstr( cntvect, lsq_2dlagrstr ) )
128 coi_error = max( coi_error, coidef_2dlagrval( cntvect, lsq_2dlagrval ) )
129
130#if defined(CONOPT_LICENSE_INT_1) && defined(CONOPT_LICENSE_INT_2) && defined(CONOPT_LICENSE_INT_3) && defined(CONOPT_LICENSE_TEXT)
131 coi_error = max( coi_error, coidef_license( cntvect, conopt_license_int_1, conopt_license_int_2, conopt_license_int_3, conopt_license_text) )
132#endif
133
134 If ( coi_error .ne. 0 ) THEN
135 write(*,*)
136 write(*,*) '**** Fatal Error while loading CONOPT Callback routines.'
137 write(*,*)
138 call flog( "Skipping Solve due to setup errors", 1 )
139 ENDIF
140!
141! Start CONOPT:
142!
143 time0 = omp_get_wtime()
144 coi_error = coi_solve( cntvect )
145 time1 = omp_get_wtime() - time0
146
147 write(*,*)
148 write(*,*) 'One thread: End of Least Square example 2. Return code=',coi_error
149
150 If ( coi_error /= 0 ) then
151 call flog( "One thread: Errors encountered during solution", 1 )
152 elseif ( stacalls == 0 .or. solcalls == 0 ) then
153 call flog( "One thread: Status or Solution routine was not called", 1 )
154 elseif ( sstat /= 1 .or. mstat /= 2 ) then
155 call flog( "One thread: Solver or Model status was not as expected (1,2)", 1 )
156 elseif ( abs( obj - 19.44434311d0 ) > 1.d-7 ) then
157 call flog( "One thread: Incorrect objective returned", 1 )
158 endif
159!
160! Start CONOPT again this time with threads
161!
162 coi_error = max( coi_error, coidef_threads( cntvect, 0 ) ) ! 0 means use as many threads as you can
163 time0 = omp_get_wtime()
164 coi_error = coi_solve( cntvect )
165 time2 = omp_get_wtime() - time0
166
167 write(*,*)
168 write(*,*) 'Many threads: End of Least Square example 2. Return code=',coi_error
169
170 If ( coi_error /= 0 ) then
171 call flog( "Many threads: Errors encountered during solution", 1 )
172 elseif ( stacalls == 0 .or. solcalls == 0 ) then
173 call flog( "Many threads: Status or Solution routine was not called", 1 )
174 elseif ( sstat /= 1 .or. mstat /= 2 ) then
175 call flog( "Many threads: Solver or Model status was not as expected (1,2)", 1 )
176 elseif ( abs( obj - 19.44434311d0 ) > 1.d-7 ) then
177 call flog( "Many threads: Incorrect objective returned", 1 )
178 endif
179
180 if ( coi_free( cntvect ) /= 0 ) call flog( "Error while freeing control vector", 1 )
181
182 write(*,*)
183 write(*,"('Time for single thread',f10.3)") time1
184 write(*,"('Time for multi thread',f10.3)") time2
185 write(*,"('Speedup ',f10.3)") time1/time2
186 write(*,"('Efficiency ',f10.3)") time1/time2/omp_get_max_threads()
187
188 call flog( "Successful Solve", 0 )
189
190End Program leastsquare
191!
192! ============================================================================
193! Define information about the model:
194!
195
196!> Define information about the model
197!!
198!! @include{doc} readMatrix_params.dox
199Integer Function lsq_readmatrix( lower, curr, upper, vsta, type, rhs, esta, &
200 colsta, rowno, value, nlflag, n, m, nz, &
201 usrmem )
202#ifdef dec_directives_win32
203!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_ReadMatrix
204#endif
205 Use lsq2
206 implicit none
207 integer, intent (in) :: n ! number of variables
208 integer, intent (in) :: m ! number of constraints
209 integer, intent (in) :: nz ! number of nonzeros
210 real*8, intent (in out), dimension(n) :: lower ! vector of lower bounds
211 real*8, intent (in out), dimension(n) :: curr ! vector of initial values
212 real*8, intent (in out), dimension(n) :: upper ! vector of upper bounds
213 integer, intent (in out), dimension(n) :: vsta ! vector of initial variable status
214 ! (not defined here)
215 integer, intent (out), dimension(m) :: type ! vector of equation types
216 integer, intent (in out), dimension(m) :: esta ! vector of initial equation status
217 ! (not defined here)
218 real*8, intent (in out), dimension(m) :: rhs ! vector of right hand sides
219 integer, intent (in out), dimension(n+1) :: colsta ! vector with start of column indices
220 integer, intent (out), dimension(nz) :: rowno ! vector of row numbers
221 integer, intent (in out), dimension(nz) :: nlflag ! vector of nonlinearity flags
222 real*8, intent (in out), dimension(nz) :: value ! vector of matrix values
223 real*8 usrmem(*) ! optional user memory
224
225 Integer :: i, j, k
226!
227! Information about Variables:
228! Default: Lower = -Inf, Curr = 0, and Upper = +inf.
229! Default: the status information in Vsta is not used.
230!
231 do i = 1, dimx
232 curr(i) = -0.8d0
233 enddo
234!
235! Information about Constraints:
236! Default: Rhs = 0
237! Default: the status information in Esta and the function
238! value in FV are not used.
239! Default: Type: There is no default.
240! 0 = Equality,
241! 1 = Greater than or equal,
242! 2 = Less than or equal,
243! 3 = Non binding.
244!
245! Constraints 1 to Nobs:
246! Rhs = Obs(i) and type Equality
247!
248 do i = 1, nobs
249 rhs(i) = obs(i)
250 type(i) = 0
251 enddo
252!
253! Constraint Nobs + 1 (Objective)
254! Rhs = 0 and type Non binding
255!
256 type(nobs+1) = 3
257!
258! Information about the Jacobian. CONOPT expects a columnwise
259! representation in Rowno, Value, Nlflag and Colsta.
260!
261! Colsta = Start of column indices (No Defaults):
262! Rowno = Row indices
263! Value = Value of derivative (by default only linear
264! derivatives are used)
265! Nlflag = 0 for linear and 1 for nonlinear derivative
266! (not needed for completely linear models)
267!
268!
269! Indices
270! x(i) res(j)
271! j: NL L=1
272! obj: NL
273!
274 k = 1
275 do j = 1, dimx
276 colsta(j) = k
277 do i = 1, nobs
278 rowno(k) = i
279 nlflag(k) = 1
280 k = k + 1
281 enddo
282 enddo
283 do i = 1, nobs
284 colsta(dimx+i) = k
285 rowno(k) = i
286 nlflag(k) = 0
287 value(k) = 1.d0
288 k = k + 1
289 rowno(k) = nobs+1
290 nlflag(k) = 1
291 k = k + 1
292 enddo
293 colsta(dimx+nobs+1) = k
294
295 lsq_readmatrix = 0 ! Return value means OK
296
297end Function lsq_readmatrix
298!
299!==========================================================================
300! Compute nonlinear terms and non-constant Jacobian elements
301!
302
303!> Compute nonlinear terms and non-constant Jacobian elements
304!!
305!! @include{doc} fdeval_params.dox
306Integer Function lsq_fdeval( x, g, jac, rowno, jcnm, mode, ignerr, errcnt, &
307 n, nz, thread, usrmem )
308#ifdef dec_directives_win32
309!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_FDEval
310#endif
311 use lsq2
312 implicit none
313 integer, intent (in) :: n ! number of variables
314 integer, intent (in) :: rowno ! number of the row to be evaluated
315 integer, intent (in) :: nz ! number of nonzeros in this row
316 real*8, intent (in), dimension(n) :: x ! vector of current solution values
317 real*8, intent (in out) :: g ! constraint value
318 real*8, intent (in out), dimension(n) :: jac ! vector of derivatives for current constraint
319 integer, intent (in), dimension(nz) :: jcnm ! list of variables that appear nonlinearly
320 ! in this row. Ffor information only.
321 integer, intent (in) :: mode ! evaluation mode: 1 = function value
322 ! 2 = derivatives, 3 = both
323 integer, intent (in) :: ignerr ! if 1 then errors can be ignored as long
324 ! as errcnt is incremented
325 integer, intent (in out) :: errcnt ! error counter to be incremented in case
326 ! of function evaluation errors.
327 integer, intent (in) :: thread
328 real*8 usrmem(*) ! optional user memory
329
330 integer :: i,j
331 real*8 :: s
332!
333! Row Nobs+1: the objective function is nonlinear
334!
335 if ( rowno .eq. nobs+1 ) then
336!
337! Mode = 1 or 3. Function value: G = sum(i, res(i)**2)
338!
339 if ( mode .eq. 1 .or. mode .eq. 3 ) then
340 s = 0.d0
341 do i = 1, nobs
342 s = s + x(dimx+i)**2
343 enddo
344 g = s
345 endif
346!
347! Mode = 2 or 3: Derivative values:
348!
349 if ( mode .eq. 2 .or. mode .eq. 3 ) then
350 do i = 1, nobs
351 jac(dimx+i) = 2*x(dimx+i)
352 enddo
353 endif
354!
355! Row 1 to Nobs: The observation definitions:
356!
357 else
358!
359! Mode = 1 or 3: Function value - x-part only
360!
361 if ( mode .eq. 1 .or. mode .eq. 3 ) then
362 s = 0.d0
363 do j = 1, dimx
364 s = s + a(rowno,j)*x(j) + b(rowno,j)*x(j)**2
365 enddo
366 g = s
367 endif
369! Mode = 2 or 3: Derivatives
370!
371 if ( mode .eq. 2 .or. mode .eq. 3 ) then
372 do j = 1, dimx
373 jac(j) = a(rowno,j) + 2.d0*b(rowno,j)*x(j)
374 enddo
375 endif
376 endif
377 lsq_fdeval = 0
378
379end Function lsq_fdeval
380
381
382!> Specify the structure of the Lagrangian of the Hessian
383!!
384!! @include{doc} 2DLagrStr_params.dox
385Integer Function lsq_2dlagrstr( HSRW, HSCL, NODRV, N, M, NHESS, UsrMem )
386#ifdef dec_directives_win32
387!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_2DLagrStr
388#endif
389 use lsq2
390 Implicit None
391
392 Integer, Intent (IN) :: n, m, nhess
393 Integer, Intent (IN OUT) :: nodrv
394 Integer, Dimension(Nhess), Intent (Out) :: hsrw, hscl
395 real*8, Intent(IN OUT) :: usrmem(*)
396
397 integer :: i
398!
399! We assume that this is a Large Residual model, where the most
400! important 2nd derivatives are those that appear directly in the
401! objective.
402! We will therefore only define the 2nd derivatives corresponding to
403! the objective constraint where the 2nd derivatives is 2*I.
404! CONOPT will by default complain that some nonlinear variables do
405! not appear in the Hessian. We will therefore define one dummy element
406! on the diagonal of the Hessian for each of the nonlinear X-variables
407! and give it the numerical value 0.
408!
409!
410! Define structure of Hessian
411!
412 do i = 1, dimx+nobs
413 hsrw(i) = i
414 hscl(i) = i
415 enddo
416
417 lsq_2dlagrstr = 0
418
419End Function lsq_2dlagrstr
420
421
422!> Compute the Lagrangian of the Hessian
423!!
424!! @include{doc} 2DLagrVal_params.dox
425Integer Function lsq_2dlagrval( X, U, HSRW, HSCL, HSVL, NODRV, N, M, NHESS, UsrMem )
426#ifdef dec_directives_win32
427!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_2DLagrVal
428#endif
429 use lsq2
430 Implicit None
431
432 Integer, Intent (IN) :: n, m, nhess
433 Integer, Intent (IN OUT) :: nodrv
434 real*8, Dimension(N), Intent (IN) :: x
435 real*8, Dimension(M), Intent (IN) :: u
436 Integer, Dimension(Nhess), Intent (In) :: hsrw, hscl
437 real*8, Dimension(NHess), Intent (Out) :: hsvl
438 real*8, Intent(IN OUT) :: usrmem(*)
439
440 integer :: i
441!
442! We assume that this is a Large Residual model, where the most
443! important 2nd derivatives are those that appear directly in the
444! objective.
445! We will therefore only define the 2nd derivatives corresponding to
446! the objective constraint where the 2nd derivatives is 2*I.
447! CONOPT will by default complain that some nonlinear variables do
448! not appear in the Hessian. We will therefore define one dummy element
449! on the diagonal of the Hessian for each of the nonlinear X-variables
450! and give it the numerical value 0.
451!
452! Normal Evaluation mode
453!
454 do i = 1, dimx
455 hsvl(i) = 0.d0
456 enddo
457 do i = 1, nobs
458 hsvl(dimx+i) = 2.d0*u(nobs+1)
459 enddo
460
461 lsq_2dlagrval = 0
462
463End Function lsq_2dlagrval
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_2dlagrstr(cntvect, coi_2dlagrstr)
define callback routine for providing the structure of the second derivatives of the Lagrangian.
Definition conopt.f90:1547
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_2dlagrval(cntvect, coi_2dlagrval)
define callback routine for computing the values of the second derivatives of the Lagrangian.
Definition conopt.f90:1573
integer(c_int) function coidef_license(cntvect, licint1, licint2, licint3, licstring)
define the License Information.
Definition conopt.f90:293
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_numhess(cntvect, numhess)
defines the Number of Hessian Nonzeros.
Definition conopt.f90:188
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_2dlagrstr(hsrw, hscl, nodrv, n, m, nhess, usrmem)
Specify the structure of the Lagrangian of the Hessian.
Definition leastsq2.f90:359
integer function lsq_2dlagrval(x, u, hsrw, hscl, hsvl, nodrv, n, m, nhess, usrmem)
Compute the Lagrangian of the Hessian.
Definition leastsq2.f90:396
#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
real *8, dimension(nobs, dimx) a
real *8, dimension(nobs, dimx) b
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