CONOPT
Loading...
Searching...
No Matches
mp_leastsq.f90
Go to the documentation of this file.
1!> @file mp_leastsq.f90
2!! @ingroup FORTOPENMP_EXAMPLES
3!!
4!! @copydoc leastsq.f90
5
6#if defined(_WIN32) && !defined(_WIN64)
7#define dec_directives_win32
8#endif
9
10module lsq
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 lsq
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 lsq
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 lsq
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#ifdef dec_directives_win32
78!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_ReadMatrix
79!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_FDEval
80!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Status
81!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Solution
82!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Message
83!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_ErrMsg
84#endif
85!
86! Control vector
87!
88 INTEGER, Dimension(:), Pointer :: cntvect
89 INTEGER :: coi_error
90
91 real*8 time0, time1, time2
92
93 call startup
94!
95! Define data
96!
97 Call defdata
98!
99! Create and initialize a Control Vector
100!
101 coi_error = coi_create( cntvect )
102!
103! Tell CONOPT about the size of the model by populating the Control Vector:
104!
105 coi_error = max( coi_error, coidef_numvar( cntvect, nobs + dimx ) )
106 coi_error = max( coi_error, coidef_numcon( cntvect, nobs + 1 ) )
107 coi_error = max( coi_error, coidef_numnz( cntvect, nobs * (dimx+2) ) )
108 coi_error = max( coi_error, coidef_numnlnz( cntvect, nobs * (dimx+1) ) )
109 coi_error = max( coi_error, coidef_optdir( cntvect, -1 ) ) ! Minimize
110 coi_error = max( coi_error, coidef_objcon( cntvect, nobs + 1 ) ) ! Objective is last constraint
111 coi_error = max( coi_error, coidef_optfile( cntvect, 'mp_leastsq.opt' ) )
112!
113! Tell CONOPT about the callback routines:
114!
115 coi_error = max( coi_error, coidef_readmatrix( cntvect, lsq_readmatrix ) )
116 coi_error = max( coi_error, coidef_fdeval( cntvect, lsq_fdeval ) )
117 coi_error = max( coi_error, coidef_status( cntvect, std_status ) )
118 coi_error = max( coi_error, coidef_solution( cntvect, std_solution ) )
119 coi_error = max( coi_error, coidef_message( cntvect, std_message ) )
120 coi_error = max( coi_error, coidef_errmsg( cntvect, std_errmsg ) )
121 coi_error = max( coi_error, coidef_debugfv( cntvect, 0 ) ) ! Debug Fdeval on or off
122
123#if defined(CONOPT_LICENSE_INT_1) && defined(CONOPT_LICENSE_INT_2) && defined(CONOPT_LICENSE_INT_3) && defined(CONOPT_LICENSE_TEXT)
124 coi_error = max( coi_error, coidef_license( cntvect, conopt_license_int_1, conopt_license_int_2, conopt_license_int_3, conopt_license_text) )
125#endif
126
127 If ( coi_error .ne. 0 ) THEN
128 write(*,*)
129 write(*,*) '**** Fatal Error while loading CONOPT Callback routines.'
130 write(*,*)
131 call flog( "Skipping Solve due to setup errors", 1 )
132 ENDIF
133!
134! Start CONOPT -- First time using a single thread:
135!
136 time0 = omp_get_wtime()
137 coi_error = coi_solve( cntvect )
138 time1 = omp_get_wtime() - time0
139
140 write(*,*)
141 write(*,*) 'End of Least Square example 1 with 1 thread. Return code=',coi_error
142
143 If ( coi_error /= 0 ) then
144 call flog( "One thread: Errors encountered during solution", 1 )
145 elseif ( stacalls == 0 .or. solcalls == 0 ) then
146 call flog( "One thread: Status or Solution routine was not called", 1 )
147 elseif ( .not. ( sstat == 1 .and. mstat == 2 ) ) then
148 call flog( "One thread: Solver or Model status was not as expected (1,2)", 1 )
149 elseif ( abs( obj - 19.44434311d0 ) > 1.d-7 ) then
150 call flog( "One thread: Incorrect objective returned", 1 )
151 endif
152!
153! Start CONOPT again -- this time using multiple threads:
154!
155 coi_error = max( coi_error, coidef_threads( cntvect, 0 ) ) ! 0 means use as many threads as you can
156 time0 = omp_get_wtime()
157 coi_error = coi_solve( cntvect )
158 time2 = omp_get_wtime() - time0
159
160 write(*,*)
161 write(*,*) 'Multi thread: End of Least Square example 1. Return code=',coi_error
162
163 If ( coi_error /= 0 ) then
164 call flog( "Multi thread: Errors encountered during solution", 1 )
165 elseif ( stacalls == 0 .or. solcalls == 0 ) then
166 call flog( "Multi thread: Status or Solution routine was not called", 1 )
167 elseif ( .not. ( sstat == 1 .and. mstat == 2 ) ) then
168 call flog( "Multi thread: Solver or Model status was not as expected (1,2)", 1 )
169 elseif ( abs( obj - 19.44434311d0 ) > 1.d-7 ) then
170 call flog( "Multi thread: Incorrect objective returned", 1 )
171 endif
172
173 if ( coi_free( cntvect ) /= 0 ) call flog( "Error while freeing control vector", 1 )
174
175 write(*,*)
176 write(*,"('Time for single thread',f10.3)") time1
177 write(*,"('Time for multi thread',f10.3)") time2
178 write(*,"('Speedup ',f10.3)") time1/time2
179 write(*,"('Efficiency ',f10.3)") time1/time2/omp_get_max_threads()
180
181 call flog( "Successful Solve", 0 )
182
183End Program leastsquare
184!
185! ============================================================================
186! Define information about the model:
187!
188
189!> Define information about the model
190!!
191!! @include{doc} readMatrix_params.dox
192Integer Function lsq_readmatrix( lower, curr, upper, vsta, type, rhs, esta, &
193 colsta, rowno, value, nlflag, n, m, nz, &
194 usrmem )
195#ifdef dec_directives_win32
196!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_ReadMatrix
197#endif
198 Use lsq
199 implicit none
200 integer, intent (in) :: n ! number of variables
201 integer, intent (in) :: m ! number of constraints
202 integer, intent (in) :: nz ! number of nonzeros
203 real*8, intent (in out), dimension(n) :: lower ! vector of lower bounds
204 real*8, intent (in out), dimension(n) :: curr ! vector of initial values
205 real*8, intent (in out), dimension(n) :: upper ! vector of upper bounds
206 integer, intent (in out), dimension(n) :: vsta ! vector of initial variable status
207 ! (not defined here)
208 integer, intent (out), dimension(m) :: type ! vector of equation types
209 integer, intent (in out), dimension(m) :: esta ! vector of initial equation status
210 ! (not defined here)
211 real*8, intent (in out), dimension(m) :: rhs ! vector of right hand sides
212 integer, intent (in out), dimension(n+1) :: colsta ! vector with start of column indices
213 integer, intent (out), dimension(nz) :: rowno ! vector of row numbers
214 integer, intent (in out), dimension(nz) :: nlflag ! vector of nonlinearity flags
215 real*8, intent (in out), dimension(nz) :: value ! vector of matrix values
216 real*8 usrmem(*) ! optional user memory
217
218 Integer :: i, j, k
219!
220! Information about Variables:
221! Default: Lower = -Inf, Curr = 0, and Upper = +inf.
222! Default: the status information in Vsta is not used.
223!
224 do i = 1, dimx
225 curr(i) = -0.8d0
226 enddo
227!
228! Information about Constraints:
229! Default: Rhs = 0
230! Default: the status information in Esta and the function
231! value in FV are not used.
232! Default: Type: There is no default.
233! 0 = Equality,
234! 1 = Greater than or equal,
235! 2 = Less than or equal,
236! 3 = Non binding.
237!
238! Constraints 1 to Nobs:
239! Rhs = Obs(i) and type Equality
240!
241 do i = 1, nobs
242 rhs(i) = obs(i)
243 type(i) = 0
244 enddo
245!
246! Constraint Nobs + 1 (Objective)
247! Rhs = 0 and type Non binding
248!
249 type(nobs+1) = 3
250!
251! Information about the Jacobian. CONOPT expects a columnwise
252! representation in Rowno, Value, Nlflag and Colsta.
253!
254! Colsta = Start of column indices (No Defaults):
255! Rowno = Row indices
256! Value = Value of derivative (by default only linear
257! derivatives are used)
258! Nlflag = 0 for linear and 1 for nonlinear derivative
259! (not needed for completely linear models)
260!
261!
262! Indices
263! x(i) res(j)
264! j: NL L=1
265! obj: NL
266!
267 k = 1
268 do j = 1, dimx
269 colsta(j) = k
270 do i = 1, nobs
271 rowno(k) = i
272 nlflag(k) = 1
273 k = k + 1
274 enddo
275 enddo
276 do i = 1, nobs
277 colsta(dimx+i) = k
278 rowno(k) = i
279 nlflag(k) = 0
280 value(k) = 1.d0
281 k = k + 1
282 rowno(k) = nobs+1
283 nlflag(k) = 1
284 k = k + 1
285 enddo
286 colsta(dimx+nobs+1) = k
288 lsq_readmatrix = 0 ! Return value means OK
289
290end Function lsq_readmatrix
291!
292!==========================================================================
293! Compute nonlinear terms and non-constant Jacobian elements
294!
295
296!> Compute nonlinear terms and non-constant Jacobian elements
297!!
298!! @include{doc} fdeval_params.dox
299Integer Function lsq_fdeval( x, g, jac, rowno, jcnm, mode, ignerr, errcnt, &
300 n, nz, thread, usrmem )
301#ifdef dec_directives_win32
302!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_FDEval
303#endif
304 use lsq
305 implicit none
306 integer, intent (in) :: n ! number of variables
307 integer, intent (in) :: rowno ! number of the row to be evaluated
308 integer, intent (in) :: nz ! number of nonzeros in this row
309 real*8, intent (in), dimension(n) :: x ! vector of current solution values
310 real*8, intent (in out) :: g ! constraint value
311 real*8, intent (in out), dimension(n) :: jac ! vector of derivatives for current constraint
312 integer, intent (in), dimension(nz) :: jcnm ! list of variables that appear nonlinearly
313 ! in this row. Ffor information only.
314 integer, intent (in) :: mode ! evaluation mode: 1 = function value
315 ! 2 = derivatives, 3 = both
316 integer, intent (in) :: ignerr ! if 1 then errors can be ignored as long
317 ! as errcnt is incremented
318 integer, intent (in out) :: errcnt ! error counter to be incremented in case
319 ! of function evaluation errors.
320 integer, intent (in) :: thread
321 real*8 usrmem(*) ! optional user memory
322
323 integer :: i,j
324 real*8 :: s
325!
326! Row Nobs+1: the objective function is nonlinear
327!
328 if ( rowno .eq. nobs+1 ) then
329!
330! Mode = 1 or 3. Function value: G = sum(i, res(i)**2)
331!
332 if ( mode .eq. 1 .or. mode .eq. 3 ) then
333 s = 0.d0
334 do i = 1, nobs
335 s = s + x(dimx+i)**2
336 enddo
337 g = s
338 endif
339!
340! Mode = 2 or 3: Derivative values:
341!
342 if ( mode .eq. 2 .or. mode .eq. 3 ) then
343 do i = 1, nobs
344 jac(dimx+i) = 2*x(dimx+i)
345 enddo
346 endif
347!
348! Row 1 to Nobs: The observation definitions:
349!
350 else
351!
352! Mode = 1 or 3: Function value - x-part only
353!
354 if ( mode .eq. 1 .or. mode .eq. 3 ) then
355 s = 0.d0
356 do j = 1, dimx
357 s = s + a(rowno,j)*x(j) + b(rowno,j)*x(j)**2
358 enddo
359 g = s
360 endif
361!
362! Mode = 2 or 3: Derivatives
363!
364 if ( mode .eq. 2 .or. mode .eq. 3 ) then
365 do j = 1, dimx
366 jac(j) = a(rowno,j) + 2.d0*b(rowno,j)*x(j)
367 enddo
368 endif
369 endif
370 lsq_fdeval = 0
371
372end Function lsq_fdeval
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_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_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
real *8, dimension(nobs) obs
real *8, dimension(nobs, dimx) b
real *8, dimension(nobs, dimx) a
subroutine flog(msg, code)
Definition comdeclp.f90:48
real *8 obj
Definition comdeclp.f90:19
subroutine startup
Definition comdeclp.f90:31